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");
320 if (!fr->bExcl_IntraCGAll_InterCGNone)
322 gmx_fatal(FARGS,"The charge-group - charge-group force loops only support systems with all intra-cg interactions excluded and no inter-cg exclusions, this is not the case for this system.");
326 if (fr->solvent_opt == esolTIP4P) {
327 enlist_w = enlistTIP4P_ATOM;
328 enlist_ww = enlistTIP4P_TIP4P;
330 enlist_w = enlistSPC_ATOM;
331 enlist_ww = enlistSPC_SPC;
334 for(i=0; i<fr->nnblists; i++)
336 nbl = &(fr->nblists[i]);
337 init_nblist(&nbl->nlist_sr[eNL_VDWQQ],&nbl->nlist_lr[eNL_VDWQQ],
338 maxsr,maxlr,ivdw,icoul,FALSE,enlist_def);
339 init_nblist(&nbl->nlist_sr[eNL_VDW],&nbl->nlist_lr[eNL_VDW],
340 maxsr,maxlr,ivdw,0,FALSE,enlist_def);
341 init_nblist(&nbl->nlist_sr[eNL_QQ],&nbl->nlist_lr[eNL_QQ],
342 maxsr,maxlr,0,icoul,FALSE,enlist_def);
343 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_WATER],&nbl->nlist_lr[eNL_VDWQQ_WATER],
344 maxsr_wat,maxlr_wat,ivdw,icoul, FALSE,enlist_w);
345 init_nblist(&nbl->nlist_sr[eNL_QQ_WATER],&nbl->nlist_lr[eNL_QQ_WATER],
346 maxsr_wat,maxlr_wat,0,icoul, FALSE,enlist_w);
347 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_WATERWATER],&nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
348 maxsr_wat,maxlr_wat,ivdw,icoul, FALSE,enlist_ww);
349 init_nblist(&nbl->nlist_sr[eNL_QQ_WATERWATER],&nbl->nlist_lr[eNL_QQ_WATERWATER],
350 maxsr_wat,maxlr_wat,0,icoul, FALSE,enlist_ww);
352 if (fr->efep != efepNO)
363 init_nblist(&nbl->nlist_sr[eNL_VDWQQ_FREE],&nbl->nlist_lr[eNL_VDWQQ_FREE],
364 maxsr,maxlr,ivdw,icoulf,TRUE,enlistATOM_ATOM);
365 init_nblist(&nbl->nlist_sr[eNL_VDW_FREE],&nbl->nlist_lr[eNL_VDW_FREE],
366 maxsr,maxlr,ivdw,0,TRUE,enlistATOM_ATOM);
367 init_nblist(&nbl->nlist_sr[eNL_QQ_FREE],&nbl->nlist_lr[eNL_QQ_FREE],
368 maxsr,maxlr,0,icoulf,TRUE,enlistATOM_ATOM);
372 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
374 init_nblist(&fr->QMMMlist,NULL,
375 maxsr,maxlr,0,icoul,FALSE,enlistATOM_ATOM);
378 fr->ns.nblist_initialized=TRUE;
381 static void reset_nblist(t_nblist *nl)
392 static void reset_neighbor_list(t_forcerec *fr,gmx_bool bLR,int nls,int eNL)
398 reset_nblist(&(fr->nblists[nls].nlist_lr[eNL]));
402 for(n=0; n<fr->nnblists; n++)
404 for(i=0; i<eNL_NR; i++)
406 reset_nblist(&(fr->nblists[n].nlist_sr[i]));
411 /* only reset the short-range nblist */
412 reset_nblist(&(fr->QMMMlist));
420 static inline void new_i_nblist(t_nblist *nlist,
421 gmx_bool bLR,atom_id i_atom,int shift,int gid)
427 /* Check whether we have to increase the i counter */
429 (nlist->iinr[nri] != i_atom) ||
430 (nlist->shift[nri] != shift) ||
431 (nlist->gid[nri] != gid))
433 /* This is something else. Now see if any entries have
434 * been added in the list of the previous atom.
437 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
438 (nlist->gid[nri] != -1)))
440 /* If so increase the counter */
443 if (nlist->nri >= nlist->maxnri)
445 nlist->maxnri += over_alloc_large(nlist->nri);
446 reallocate_nblist(nlist);
449 /* Set the number of neighbours and the atom number */
450 nlist->jindex[nri+1] = nlist->jindex[nri];
451 nlist->iinr[nri] = i_atom;
452 nlist->gid[nri] = gid;
453 nlist->shift[nri] = shift;
457 static inline void close_i_nblist(t_nblist *nlist)
459 int nri = nlist->nri;
464 nlist->jindex[nri+1] = nlist->nrj;
466 len=nlist->nrj - nlist->jindex[nri];
468 /* nlist length for water i molecules is treated statically
471 if (len > nlist->maxlen)
478 static inline void close_nblist(t_nblist *nlist)
480 /* Only close this nblist when it has been initialized */
487 static inline void close_neighbor_list(t_forcerec *fr,gmx_bool bLR,int nls,int eNL,
488 gmx_bool bMakeQMMMnblist)
492 if (bMakeQMMMnblist) {
495 close_nblist(&(fr->QMMMlist));
502 close_nblist(&(fr->nblists[nls].nlist_lr[eNL]));
506 for(n=0; n<fr->nnblists; n++)
508 for(i=0; (i<eNL_NR); i++)
510 close_nblist(&(fr->nblists[n].nlist_sr[i]));
517 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
521 if (nlist->nrj >= nlist->maxnrj)
523 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
525 fprintf(debug,"Increasing %s nblist %s j size to %d\n",
526 bLR ? "LR" : "SR",nrnb_str(nlist->il_code),nlist->maxnrj);
528 srenew(nlist->jjnr,nlist->maxnrj);
531 nlist->jjnr[nrj] = j_atom;
535 static inline void add_j_to_nblist_cg(t_nblist *nlist,
536 atom_id j_start,int j_end,
537 t_excl *bexcl,gmx_bool bLR)
542 if (nlist->nrj >= nlist->maxnrj)
544 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
546 fprintf(debug,"Increasing %s nblist %s j size to %d\n",
547 bLR ? "LR" : "SR",nrnb_str(nlist->il_code),nlist->maxnrj);
549 srenew(nlist->jjnr ,nlist->maxnrj);
550 srenew(nlist->jjnr_end,nlist->maxnrj);
551 srenew(nlist->excl ,nlist->maxnrj*MAX_CGCGSIZE);
554 nlist->jjnr[nrj] = j_start;
555 nlist->jjnr_end[nrj] = j_end;
557 if (j_end - j_start > MAX_CGCGSIZE)
559 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);
562 /* Set the exclusions */
563 for(j=j_start; j<j_end; j++)
565 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
572 put_in_list_t(gmx_bool bHaveVdW[],
588 put_in_list_at(gmx_bool bHaveVdW[],
603 /* The a[] index has been removed,
604 * to put it back in i_atom should be a[i0] and jj should be a[jj].
609 t_nblist * vdwc_free = NULL;
610 t_nblist * vdw_free = NULL;
611 t_nblist * coul_free = NULL;
612 t_nblist * vdwc_ww = NULL;
613 t_nblist * coul_ww = NULL;
615 int i,j,jcg,igid,gid,nbl_ind,ind_ij;
616 atom_id jj,jj0,jj1,i_atom;
621 real *charge,*chargeB;
623 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
624 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
628 /* Copy some pointers */
630 charge = md->chargeA;
631 chargeB = md->chargeB;
634 bPert = md->bPerturbed;
638 nicg = index[icg+1]-i0;
640 /* Get the i charge group info */
641 igid = GET_CGINFO_GID(cginfo[icg]);
642 iwater = GET_CGINFO_SOLOPT(cginfo[icg]);
647 /* Check if any of the particles involved are perturbed.
648 * If not we can do the cheaper normal put_in_list
649 * and use more solvent optimization.
651 for(i=0; i<nicg; i++)
653 bFreeEnergy |= bPert[i0+i];
655 /* Loop over the j charge groups */
656 for(j=0; (j<nj && !bFreeEnergy); j++)
661 /* Finally loop over the atoms in the j-charge group */
662 for(jj=jj0; jj<jj1; jj++)
664 bFreeEnergy |= bPert[jj];
669 /* Unpack pointers to neighbourlist structs */
670 if (fr->nnblists == 1)
676 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
680 nlist = fr->nblists[nbl_ind].nlist_lr;
684 nlist = fr->nblists[nbl_ind].nlist_sr;
687 if (iwater != esolNO)
689 vdwc = &nlist[eNL_VDWQQ_WATER];
690 vdw = &nlist[eNL_VDW];
691 coul = &nlist[eNL_QQ_WATER];
692 #ifndef DISABLE_WATERWATER_NLIST
693 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
694 coul_ww = &nlist[eNL_QQ_WATERWATER];
699 vdwc = &nlist[eNL_VDWQQ];
700 vdw = &nlist[eNL_VDW];
701 coul = &nlist[eNL_QQ];
706 if (iwater != esolNO)
708 /* Loop over the atoms in the i charge group */
710 gid = GID(igid,jgid,ngid);
711 /* Create new i_atom for each energy group */
712 if (bDoCoul && bDoVdW)
714 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
715 #ifndef DISABLE_WATERWATER_NLIST
716 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
721 new_i_nblist(vdw,bLR,i_atom,shift,gid);
725 new_i_nblist(coul,bLR,i_atom,shift,gid);
726 #ifndef DISABLE_WATERWATER_NLIST
727 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
730 /* Loop over the j charge groups */
731 for(j=0; (j<nj); j++)
741 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
743 if (iwater == esolSPC && jwater == esolSPC)
745 /* Interaction between two SPC molecules */
748 /* VdW only - only first atoms in each water interact */
749 add_j_to_nblist(vdw,jj0,bLR);
753 #ifdef DISABLE_WATERWATER_NLIST
754 /* Add entries for the three atoms - only do VdW if we need to */
757 add_j_to_nblist(coul,jj0,bLR);
761 add_j_to_nblist(vdwc,jj0,bLR);
763 add_j_to_nblist(coul,jj0+1,bLR);
764 add_j_to_nblist(coul,jj0+2,bLR);
766 /* One entry for the entire water-water interaction */
769 add_j_to_nblist(coul_ww,jj0,bLR);
773 add_j_to_nblist(vdwc_ww,jj0,bLR);
778 else if (iwater == esolTIP4P && jwater == esolTIP4P)
780 /* Interaction between two TIP4p molecules */
783 /* VdW only - only first atoms in each water interact */
784 add_j_to_nblist(vdw,jj0,bLR);
788 #ifdef DISABLE_WATERWATER_NLIST
789 /* Add entries for the four atoms - only do VdW if we need to */
792 add_j_to_nblist(vdw,jj0,bLR);
794 add_j_to_nblist(coul,jj0+1,bLR);
795 add_j_to_nblist(coul,jj0+2,bLR);
796 add_j_to_nblist(coul,jj0+3,bLR);
798 /* One entry for the entire water-water interaction */
801 add_j_to_nblist(coul_ww,jj0,bLR);
805 add_j_to_nblist(vdwc_ww,jj0,bLR);
812 /* j charge group is not water, but i is.
813 * Add entries to the water-other_atom lists; the geometry of the water
814 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
815 * so we don't care if it is SPC or TIP4P...
822 for(jj=jj0; (jj<jj1); jj++)
826 add_j_to_nblist(coul,jj,bLR);
832 for(jj=jj0; (jj<jj1); jj++)
834 if (bHaveVdW[type[jj]])
836 add_j_to_nblist(vdw,jj,bLR);
842 /* _charge_ _groups_ interact with both coulomb and LJ */
843 /* Check which atoms we should add to the lists! */
844 for(jj=jj0; (jj<jj1); jj++)
846 if (bHaveVdW[type[jj]])
850 add_j_to_nblist(vdwc,jj,bLR);
854 add_j_to_nblist(vdw,jj,bLR);
857 else if (charge[jj] != 0)
859 add_j_to_nblist(coul,jj,bLR);
866 close_i_nblist(coul);
867 close_i_nblist(vdwc);
868 #ifndef DISABLE_WATERWATER_NLIST
869 close_i_nblist(coul_ww);
870 close_i_nblist(vdwc_ww);
875 /* no solvent as i charge group */
876 /* Loop over the atoms in the i charge group */
877 for(i=0; i<nicg; i++)
880 gid = GID(igid,jgid,ngid);
883 /* Create new i_atom for each energy group */
884 if (bDoVdW && bDoCoul)
886 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
890 new_i_nblist(vdw,bLR,i_atom,shift,gid);
894 new_i_nblist(coul,bLR,i_atom,shift,gid);
896 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
897 bDoCoul_i = (bDoCoul && qi!=0);
899 if (bDoVdW_i || bDoCoul_i)
901 /* Loop over the j charge groups */
902 for(j=0; (j<nj); j++)
906 /* Check for large charge groups */
917 /* Finally loop over the atoms in the j-charge group */
918 for(jj=jj0; jj<jj1; jj++)
920 bNotEx = NOTEXCL(bExcl,i,jj);
928 add_j_to_nblist(coul,jj,bLR);
933 if (bHaveVdW[type[jj]])
935 add_j_to_nblist(vdw,jj,bLR);
940 if (bHaveVdW[type[jj]])
944 add_j_to_nblist(vdwc,jj,bLR);
948 add_j_to_nblist(vdw,jj,bLR);
951 else if (charge[jj] != 0)
953 add_j_to_nblist(coul,jj,bLR);
961 close_i_nblist(coul);
962 close_i_nblist(vdwc);
968 /* we are doing free energy */
969 vdwc_free = &nlist[eNL_VDWQQ_FREE];
970 vdw_free = &nlist[eNL_VDW_FREE];
971 coul_free = &nlist[eNL_QQ_FREE];
972 /* Loop over the atoms in the i charge group */
973 for(i=0; i<nicg; i++)
976 gid = GID(igid,jgid,ngid);
978 qiB = chargeB[i_atom];
980 /* Create new i_atom for each energy group */
981 if (bDoVdW && bDoCoul)
982 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
984 new_i_nblist(vdw,bLR,i_atom,shift,gid);
986 new_i_nblist(coul,bLR,i_atom,shift,gid);
988 new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
989 new_i_nblist(coul_free,bLR,i_atom,shift,gid);
990 new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
992 bDoVdW_i = (bDoVdW &&
993 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
994 bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
995 /* For TIP4P the first atom does not have a charge,
996 * but the last three do. So we should still put an atom
997 * without LJ but with charge in the water-atom neighborlist
998 * for a TIP4p i charge group.
999 * For SPC type water the first atom has LJ and charge,
1000 * so there is no such problem.
1002 if (iwater == esolNO)
1004 bDoCoul_i_sol = bDoCoul_i;
1008 bDoCoul_i_sol = bDoCoul;
1011 if (bDoVdW_i || bDoCoul_i_sol)
1013 /* Loop over the j charge groups */
1014 for(j=0; (j<nj); j++)
1018 /* Check for large charge groups */
1029 /* Finally loop over the atoms in the j-charge group */
1030 bFree = bPert[i_atom];
1031 for(jj=jj0; (jj<jj1); jj++)
1033 bFreeJ = bFree || bPert[jj];
1034 /* Complicated if, because the water H's should also
1035 * see perturbed j-particles
1037 if (iwater==esolNO || i==0 || bFreeJ)
1039 bNotEx = NOTEXCL(bExcl,i,jj);
1047 if (charge[jj]!=0 || chargeB[jj]!=0)
1049 add_j_to_nblist(coul_free,jj,bLR);
1052 else if (!bDoCoul_i)
1054 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1056 add_j_to_nblist(vdw_free,jj,bLR);
1061 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1063 if (charge[jj]!=0 || chargeB[jj]!=0)
1065 add_j_to_nblist(vdwc_free,jj,bLR);
1069 add_j_to_nblist(vdw_free,jj,bLR);
1072 else if (charge[jj]!=0 || chargeB[jj]!=0)
1073 add_j_to_nblist(coul_free,jj,bLR);
1078 /* This is done whether or not bWater is set */
1079 if (charge[jj] != 0)
1081 add_j_to_nblist(coul,jj,bLR);
1084 else if (!bDoCoul_i_sol)
1086 if (bHaveVdW[type[jj]])
1088 add_j_to_nblist(vdw,jj,bLR);
1093 if (bHaveVdW[type[jj]])
1095 if (charge[jj] != 0)
1097 add_j_to_nblist(vdwc,jj,bLR);
1101 add_j_to_nblist(vdw,jj,bLR);
1104 else if (charge[jj] != 0)
1106 add_j_to_nblist(coul,jj,bLR);
1114 close_i_nblist(vdw);
1115 close_i_nblist(coul);
1116 close_i_nblist(vdwc);
1117 close_i_nblist(vdw_free);
1118 close_i_nblist(coul_free);
1119 close_i_nblist(vdwc_free);
1125 put_in_list_qmmm(gmx_bool bHaveVdW[],
1141 int i,j,jcg,igid,gid;
1142 atom_id jj,jj0,jj1,i_atom;
1146 /* Get atom range */
1148 nicg = index[icg+1]-i0;
1150 /* Get the i charge group info */
1151 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1153 coul = &fr->QMMMlist;
1155 /* Loop over atoms in the ith charge group */
1156 for (i=0;i<nicg;i++)
1159 gid = GID(igid,jgid,ngid);
1160 /* Create new i_atom for each energy group */
1161 new_i_nblist(coul,bLR,i_atom,shift,gid);
1163 /* Loop over the j charge groups */
1168 /* Charge groups cannot have QM and MM atoms simultaneously */
1173 /* Finally loop over the atoms in the j-charge group */
1174 for(jj=jj0; jj<jj1; jj++)
1176 bNotEx = NOTEXCL(bExcl,i,jj);
1178 add_j_to_nblist(coul,jj,bLR);
1182 close_i_nblist(coul);
1187 put_in_list_cg(gmx_bool bHaveVdW[],
1203 int igid,gid,nbl_ind;
1207 cginfo = fr->cginfo[icg];
1209 igid = GET_CGINFO_GID(cginfo);
1210 gid = GID(igid,jgid,ngid);
1212 /* Unpack pointers to neighbourlist structs */
1213 if (fr->nnblists == 1)
1219 nbl_ind = fr->gid2nblists[gid];
1223 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1227 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1230 /* Make a new neighbor list for charge group icg.
1231 * Currently simply one neighbor list is made with LJ and Coulomb.
1232 * If required, zero interactions could be removed here
1233 * or in the force loop.
1235 new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1236 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1238 for(j=0; (j<nj); j++)
1241 /* Skip the icg-icg pairs if all self interactions are excluded */
1242 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1244 /* Here we add the j charge group jcg to the list,
1245 * exclusions are also added to the list.
1247 add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,bLR);
1251 close_i_nblist(vdwc);
1254 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1261 for(i=start; i<end; i++)
1263 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1265 SETEXCL(bexcl,i-start,excl->a[k]);
1271 for(i=start; i<end; i++)
1273 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1275 RMEXCL(bexcl,i-start,excl->a[k]);
1281 int calc_naaj(int icg,int cgtot)
1285 if ((cgtot % 2) == 1)
1287 /* Odd number of charge groups, easy */
1288 naaj = 1 + (cgtot/2);
1290 else if ((cgtot % 4) == 0)
1292 /* Multiple of four is hard */
1329 fprintf(log,"naaj=%d\n",naaj);
1335 /************************************************
1337 * S I M P L E C O R E S T U F F
1339 ************************************************/
1341 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1342 rvec b_inv,int *shift)
1344 /* This code assumes that the cut-off is smaller than
1345 * a half times the smallest diagonal element of the box.
1352 /* Compute diff vector */
1353 dz = xj[ZZ] - xi[ZZ];
1354 dy = xj[YY] - xi[YY];
1355 dx = xj[XX] - xi[XX];
1357 /* Perform NINT operation, using trunc operation, therefore
1358 * we first add 2.5 then subtract 2 again
1360 tz = dz*b_inv[ZZ] + h25;
1362 dz -= tz*box[ZZ][ZZ];
1363 dy -= tz*box[ZZ][YY];
1364 dx -= tz*box[ZZ][XX];
1366 ty = dy*b_inv[YY] + h25;
1368 dy -= ty*box[YY][YY];
1369 dx -= ty*box[YY][XX];
1371 tx = dx*b_inv[XX]+h25;
1373 dx -= tx*box[XX][XX];
1375 /* Distance squared */
1376 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1378 *shift = XYZ2IS(tx,ty,tz);
1383 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1384 rvec b_inv,int *shift)
1392 /* Compute diff vector */
1393 dx = xj[XX] - xi[XX];
1394 dy = xj[YY] - xi[YY];
1395 dz = xj[ZZ] - xi[ZZ];
1397 /* Perform NINT operation, using trunc operation, therefore
1398 * we first add 1.5 then subtract 1 again
1400 tx = dx*b_inv[XX] + h15;
1401 ty = dy*b_inv[YY] + h15;
1402 tz = dz*b_inv[ZZ] + h15;
1407 /* Correct diff vector for translation */
1408 ddx = tx*box_size[XX] - dx;
1409 ddy = ty*box_size[YY] - dy;
1410 ddz = tz*box_size[ZZ] - dz;
1412 /* Distance squared */
1413 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1415 *shift = XYZ2IS(tx,ty,tz);
1420 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1421 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1422 int icg,int jgid,t_block *cgs,t_excl bexcl[],
1423 int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1425 if (nsbuf->nj + nrj > MAX_CG)
1427 put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1428 cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE);
1429 /* Reset buffer contents */
1430 nsbuf->ncg = nsbuf->nj = 0;
1432 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1436 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1437 int njcg,atom_id jcg[],
1438 matrix box,rvec b_inv,real rcut2,
1439 t_block *cgs,t_ns_buf **ns_buf,
1440 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1441 t_excl bexcl[],t_forcerec *fr,
1442 put_in_list_t *put_in_list)
1446 int *cginfo=fr->cginfo;
1447 atom_id cg_j,*cgindex;
1450 cgindex = cgs->index;
1452 for(j=0; (j<njcg); j++)
1455 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1456 if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1458 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1459 if (!(i_egp_flags[jgid] & EGP_EXCL))
1461 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1462 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1469 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1470 int njcg,atom_id jcg[],
1471 gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1472 t_block *cgs,t_ns_buf **ns_buf,
1473 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1474 t_excl bexcl[],t_forcerec *fr,
1475 put_in_list_t *put_in_list)
1479 int *cginfo=fr->cginfo;
1480 atom_id cg_j,*cgindex;
1483 cgindex = cgs->index;
1487 for(j=0; (j<njcg); j++)
1490 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1491 if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1493 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1494 if (!(i_egp_flags[jgid] & EGP_EXCL))
1496 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1497 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1505 for(j=0; (j<njcg); j++)
1508 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1509 if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1510 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1511 if (!(i_egp_flags[jgid] & EGP_EXCL))
1513 add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1514 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1522 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1524 static int ns_simple_core(t_forcerec *fr,
1525 gmx_localtop_t *top,
1527 matrix box,rvec box_size,
1528 t_excl bexcl[],atom_id *aaj,
1529 int ngid,t_ns_buf **ns_buf,
1530 put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1534 int nsearch,icg,jcg,igid,i0,nri,nn;
1537 /* atom_id *i_atoms; */
1538 t_block *cgs=&(top->cgs);
1539 t_blocka *excl=&(top->excls);
1542 gmx_bool bBox,bTriclinic;
1545 rlist2 = sqr(fr->rlist);
1547 bBox = (fr->ePBC != epbcNONE);
1550 for(m=0; (m<DIM); m++)
1552 b_inv[m] = divide_err(1.0,box_size[m]);
1554 bTriclinic = TRICLINIC(box);
1561 cginfo = fr->cginfo;
1564 for (icg=fr->cg0; (icg<fr->hcg); icg++)
1567 i0 = cgs->index[icg];
1568 nri = cgs->index[icg+1]-i0;
1569 i_atoms = &(cgs->a[i0]);
1570 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1571 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1573 igid = GET_CGINFO_GID(cginfo[icg]);
1574 i_egp_flags = fr->egp_flags + ngid*igid;
1575 setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1577 naaj=calc_naaj(icg,cgs->nr);
1580 ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1581 box,b_inv,rlist2,cgs,ns_buf,
1582 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1586 ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1587 bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1588 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1592 for(nn=0; (nn<ngid); nn++)
1594 for(k=0; (k<SHIFTS); k++)
1596 nsbuf = &(ns_buf[nn][k]);
1599 put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1600 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE);
1601 nsbuf->ncg=nsbuf->nj=0;
1605 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1606 setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1608 close_neighbor_list(fr,FALSE,-1,-1,FALSE);
1613 /************************************************
1615 * N S 5 G R I D S T U F F
1617 ************************************************/
1619 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1620 int *dx0,int *dx1,real *dcx2)
1648 for(i=xgi0; i>=0; i--)
1650 dcx = (i+1)*gridx-x;
1657 for(i=xgi1; i<Nx; i++)
1670 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1671 int ncpddc,int shift_min,int shift_max,
1672 int *g0,int *g1,real *dcx2)
1675 int g_min,g_max,shift_home;
1708 g_min = (shift_min == shift_home ? 0 : ncpddc);
1709 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1716 else if (shift_max < 0)
1731 /* Check one grid cell down */
1732 dcx = ((*g0 - 1) + 1)*gridx - x;
1744 /* Check one grid cell up */
1745 dcx = (*g1 + 1)*gridx - x;
1757 #define sqr(x) ((x)*(x))
1758 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1759 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1760 /****************************************************
1762 * F A S T N E I G H B O R S E A R C H I N G
1764 * Optimized neighboursearching routine using grid
1765 * at least 1x1x1, see GROMACS manual
1767 ****************************************************/
1769 static void do_longrange(t_commrec *cr,gmx_localtop_t *top,t_forcerec *fr,
1770 int ngid,t_mdatoms *md,int icg,
1772 atom_id lr[],t_excl bexcl[],int shift,
1773 rvec x[],rvec box_size,t_nrnb *nrnb,
1774 real lambda,real *dvdlambda,
1775 gmx_grppairener_t *grppener,
1776 gmx_bool bDoVdW,gmx_bool bDoCoul,
1777 gmx_bool bEvaluateNow,put_in_list_t *put_in_list,
1778 gmx_bool bHaveVdW[],
1779 gmx_bool bDoForces,rvec *f)
1784 for(n=0; n<fr->nnblists; n++)
1786 for(i=0; (i<eNL_NR); i++)
1788 nl = &fr->nblists[n].nlist_lr[i];
1789 if ((nl->nri > nl->maxnri-32) || bEvaluateNow)
1791 close_neighbor_list(fr,TRUE,n,i,FALSE);
1792 /* Evaluate the energies and forces */
1793 do_nonbonded(cr,fr,x,f,md,NULL,
1794 grppener->ener[fr->bBHAM ? egBHAMLR : egLJLR],
1795 grppener->ener[egCOULLR],
1796 grppener->ener[egGB],box_size,
1797 nrnb,lambda,dvdlambda,n,i,
1798 GMX_DONB_LR | GMX_DONB_FORCES);
1800 reset_neighbor_list(fr,TRUE,n,i);
1807 /* Put the long range particles in a list */
1808 /* do_longrange is never called for QMMM */
1809 put_in_list(bHaveVdW,ngid,md,icg,jgid,nlr,lr,top->cgs.index,
1810 bexcl,shift,fr,TRUE,bDoVdW,bDoCoul);
1814 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1815 real *rvdw2,real *rcoul2,
1816 real *rs2,real *rm2,real *rl2)
1818 *rs2 = sqr(fr->rlist);
1819 if (bDoLongRange && fr->bTwinRange)
1821 /* The VdW and elec. LR cut-off's could be different,
1822 * so we can not simply set them to rlistlong.
1824 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1825 fr->rvdw > fr->rlist)
1827 *rvdw2 = sqr(fr->rlistlong);
1831 *rvdw2 = sqr(fr->rvdw);
1833 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1834 fr->rcoulomb > fr->rlist)
1836 *rcoul2 = sqr(fr->rlistlong);
1840 *rcoul2 = sqr(fr->rcoulomb);
1845 /* Workaround for a gcc -O3 or -ffast-math problem */
1849 *rm2 = min(*rvdw2,*rcoul2);
1850 *rl2 = max(*rvdw2,*rcoul2);
1853 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
1855 real rvdw2,rcoul2,rs2,rm2,rl2;
1858 get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1860 /* Short range buffers */
1861 snew(ns->nl_sr,ngid);
1864 snew(ns->nlr_ljc,ngid);
1865 snew(ns->nlr_one,ngid);
1869 /* Long range VdW and Coul buffers */
1870 snew(ns->nl_lr_ljc,ngid);
1874 /* Long range VdW or Coul only buffers */
1875 snew(ns->nl_lr_one,ngid);
1877 for(j=0; (j<ngid); j++) {
1878 snew(ns->nl_sr[j],MAX_CG);
1881 snew(ns->nl_lr_ljc[j],MAX_CG);
1885 snew(ns->nl_lr_one[j],MAX_CG);
1891 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
1896 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
1897 matrix box,rvec box_size,int ngid,
1898 gmx_localtop_t *top,
1899 t_grid *grid,rvec x[],
1900 t_excl bexcl[],gmx_bool *bExcludeAlleg,
1901 t_nrnb *nrnb,t_mdatoms *md,
1902 real lambda,real *dvdlambda,
1903 gmx_grppairener_t *grppener,
1904 put_in_list_t *put_in_list,
1905 gmx_bool bHaveVdW[],
1906 gmx_bool bDoLongRange,gmx_bool bDoForces,rvec *f,
1907 gmx_bool bMakeQMMMnblist)
1910 atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
1911 int *nlr_ljc,*nlr_one,*nsr;
1912 gmx_domdec_t *dd=NULL;
1913 t_block *cgs=&(top->cgs);
1914 int *cginfo=fr->cginfo;
1915 /* atom_id *i_atoms,*cgsindex=cgs->index; */
1917 int cell_x,cell_y,cell_z;
1918 int d,tx,ty,tz,dx,dy,dz,cj;
1919 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1920 int zsh_ty,zsh_tx,ysh_tx;
1922 int dx0,dx1,dy0,dy1,dz0,dz1;
1923 int Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
1924 real gridx,gridy,gridz,grid_x,grid_y,grid_z;
1925 real *dcx2,*dcy2,*dcz2;
1927 int cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
1928 int jcg0,jcg1,jjcg,cgj0,jgid;
1929 int *grida,*gridnra,*gridind;
1930 gmx_bool rvdw_lt_rcoul,rcoul_lt_rvdw;
1931 rvec xi,*cgcm,grid_offset;
1932 real r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
1934 gmx_bool bDomDec,bTriclinicX,bTriclinicY;
1939 bDomDec = DOMAINDECOMP(cr);
1945 bTriclinicX = ((YY < grid->npbcdim &&
1946 (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
1947 (ZZ < grid->npbcdim &&
1948 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
1949 bTriclinicY = (ZZ < grid->npbcdim &&
1950 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
1954 get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1956 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
1957 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
1959 if (bMakeQMMMnblist)
1967 nl_lr_ljc = ns->nl_lr_ljc;
1968 nl_lr_one = ns->nl_lr_one;
1969 nlr_ljc = ns->nlr_ljc;
1970 nlr_one = ns->nlr_one;
1978 gridind = grid->index;
1979 gridnra = grid->nra;
1982 gridx = grid->cell_size[XX];
1983 gridy = grid->cell_size[YY];
1984 gridz = grid->cell_size[ZZ];
1988 copy_rvec(grid->cell_offset,grid_offset);
1989 copy_ivec(grid->ncpddc,ncpddc);
1994 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1995 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
1996 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
1997 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
1998 if (zsh_tx!=0 && ysh_tx!=0)
2000 /* This could happen due to rounding, when both ratios are 0.5 */
2009 /* We only want a list for the test particle */
2018 /* Set the shift range */
2019 for(d=0; d<DIM; d++)
2023 /* Check if we need periodicity shifts.
2024 * Without PBC or with domain decomposition we don't need them.
2026 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2033 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2044 /* Loop over charge groups */
2045 for(icg=cg0; (icg < cg1); icg++)
2047 igid = GET_CGINFO_GID(cginfo[icg]);
2048 /* Skip this charge group if all energy groups are excluded! */
2049 if (bExcludeAlleg[igid])
2054 i0 = cgs->index[icg];
2056 if (bMakeQMMMnblist)
2058 /* Skip this charge group if it is not a QM atom while making a
2059 * QM/MM neighbourlist
2061 if (md->bQM[i0]==FALSE)
2063 continue; /* MM particle, go to next particle */
2066 /* Compute the number of charge groups that fall within the control
2069 naaj = calc_naaj(icg,cgsnr);
2076 /* make a normal neighbourlist */
2080 /* Get the j charge-group and dd cell shift ranges */
2081 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
2086 /* Compute the number of charge groups that fall within the control
2089 naaj = calc_naaj(icg,cgsnr);
2095 /* The i-particle is awlways the test particle,
2096 * so we want all j-particles
2098 max_jcg = cgsnr - 1;
2102 max_jcg = jcg1 - cgsnr;
2107 i_egp_flags = fr->egp_flags + igid*ngid;
2109 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2110 setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
2112 ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2114 /* Changed iicg to icg, DvdS 990115
2115 * (but see consistency check above, DvdS 990330)
2118 fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2119 icg,naaj,cell_x,cell_y,cell_z);
2121 /* Loop over shift vectors in three dimensions */
2122 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2124 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2125 /* Calculate range of cells in Z direction that have the shift tz */
2126 zgi = cell_z + tz*Nz;
2129 get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2131 get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2132 ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2138 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2140 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2141 /* Calculate range of cells in Y direction that have the shift ty */
2144 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2148 ygi = cell_y + ty*Ny;
2151 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2153 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2154 ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2160 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2162 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2163 /* Calculate range of cells in X direction that have the shift tx */
2166 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2170 xgi = cell_x + tx*Nx;
2173 get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2175 get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2176 ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2182 /* Get shift vector */
2183 shift=XYZ2IS(tx,ty,tz);
2185 range_check(shift,0,SHIFTS);
2187 for(nn=0; (nn<ngid); nn++)
2194 fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2195 shift,dx0,dx1,dy0,dy1,dz0,dz1);
2196 fprintf(log,"cgcm: %8.3f %8.3f %8.3f\n",cgcm[icg][XX],
2197 cgcm[icg][YY],cgcm[icg][ZZ]);
2198 fprintf(log,"xi: %8.3f %8.3f %8.3f\n",XI,YI,ZI);
2200 for (dx=dx0; (dx<=dx1); dx++)
2202 tmp1 = rl2 - dcx2[dx];
2203 for (dy=dy0; (dy<=dy1); dy++)
2205 tmp2 = tmp1 - dcy2[dy];
2208 for (dz=dz0; (dz<=dz1); dz++) {
2209 if (tmp2 > dcz2[dz]) {
2210 /* Find grid-cell cj in which possible neighbours are */
2211 cj = xyz2ci(Ny,Nz,dx,dy,dz);
2213 /* Check out how many cgs (nrj) there in this cell */
2216 /* Find the offset in the cg list */
2219 /* Check if all j's are out of range so we
2220 * can skip the whole cell.
2221 * Should save some time, especially with DD.
2224 (grida[cgj0] >= max_jcg &&
2225 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2231 for (j=0; (j<nrj); j++)
2233 jjcg = grida[cgj0+j];
2235 /* check whether this guy is in range! */
2236 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2239 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2241 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2242 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2243 /* check energy group exclusions */
2244 if (!(i_egp_flags[jgid] & EGP_EXCL))
2248 if (nsr[jgid] >= MAX_CG)
2250 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2251 nsr[jgid],nl_sr[jgid],
2252 cgs->index,/* cgsatoms, */ bexcl,
2253 shift,fr,FALSE,TRUE,TRUE);
2256 nl_sr[jgid][nsr[jgid]++]=jjcg;
2260 if (nlr_ljc[jgid] >= MAX_CG)
2262 do_longrange(cr,top,fr,ngid,md,icg,jgid,
2264 nl_lr_ljc[jgid],bexcl,shift,x,
2274 nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2278 if (nlr_one[jgid] >= MAX_CG) {
2279 do_longrange(cr,top,fr,ngid,md,icg,jgid,
2281 nl_lr_one[jgid],bexcl,shift,x,
2285 rvdw_lt_rcoul,rcoul_lt_rvdw,FALSE,
2291 nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2303 /* CHECK whether there is anything left in the buffers */
2304 for(nn=0; (nn<ngid); nn++)
2308 put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2309 cgs->index, /* cgsatoms, */ bexcl,
2310 shift,fr,FALSE,TRUE,TRUE);
2313 if (nlr_ljc[nn] > 0)
2315 do_longrange(cr,top,fr,ngid,md,icg,nn,nlr_ljc[nn],
2316 nl_lr_ljc[nn],bexcl,shift,x,box_size,nrnb,
2317 lambda,dvdlambda,grppener,TRUE,TRUE,FALSE,
2318 put_in_list,bHaveVdW,bDoForces,f);
2321 if (nlr_one[nn] > 0)
2323 do_longrange(cr,top,fr,ngid,md,icg,nn,nlr_one[nn],
2324 nl_lr_one[nn],bexcl,shift,x,box_size,nrnb,
2325 lambda,dvdlambda,grppener,
2326 rvdw_lt_rcoul,rcoul_lt_rvdw,FALSE,
2327 put_in_list,bHaveVdW,bDoForces,f);
2333 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2334 setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2336 /* Perform any left over force calculations */
2337 for (nn=0; (nn<ngid); nn++)
2341 do_longrange(cr,top,fr,0,md,icg,nn,nlr_ljc[nn],
2342 nl_lr_ljc[nn],bexcl,shift,x,box_size,nrnb,
2343 lambda,dvdlambda,grppener,
2344 TRUE,TRUE,TRUE,put_in_list,bHaveVdW,bDoForces,f);
2347 do_longrange(cr,top,fr,0,md,icg,nn,nlr_one[nn],
2348 nl_lr_one[nn],bexcl,shift,x,box_size,nrnb,
2349 lambda,dvdlambda,grppener,
2350 rvdw_lt_rcoul,rcoul_lt_rvdw,
2351 TRUE,put_in_list,bHaveVdW,bDoForces,f);
2356 /* Close off short range neighbourlists */
2357 close_neighbor_list(fr,FALSE,-1,-1,bMakeQMMMnblist);
2362 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2366 if (natoms > ns->nra_alloc)
2368 ns->nra_alloc = over_alloc_dd(natoms);
2369 srenew(ns->bexcl,ns->nra_alloc);
2370 for(i=0; i<ns->nra_alloc; i++)
2377 void init_ns(FILE *fplog,const t_commrec *cr,
2378 gmx_ns_t *ns,t_forcerec *fr,
2379 const gmx_mtop_t *mtop,
2382 int mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2386 /* Compute largest charge groups size (# atoms) */
2388 for(mt=0; mt<mtop->nmoltype; mt++) {
2389 cgs = &mtop->moltype[mt].cgs;
2390 for (icg=0; (icg < cgs->nr); icg++)
2392 nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2396 /* Verify whether largest charge group is <= max cg.
2397 * This is determined by the type of the local exclusion type
2398 * Exclusions are stored in bits. (If the type is not large
2399 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2401 maxcg = sizeof(t_excl)*8;
2402 if (nr_in_cg > maxcg)
2404 gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2408 ngid = mtop->groups.grps[egcENER].nr;
2409 snew(ns->bExcludeAlleg,ngid);
2410 for(i=0; i<ngid; i++) {
2411 ns->bExcludeAlleg[i] = TRUE;
2412 for(j=0; j<ngid; j++)
2414 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2416 ns->bExcludeAlleg[i] = FALSE;
2423 ns->grid = init_grid(fplog,fr);
2424 init_nsgrid_lists(fr,ngid,ns);
2429 snew(ns->ns_buf,ngid);
2430 for(i=0; (i<ngid); i++)
2432 snew(ns->ns_buf[i],SHIFTS);
2434 ncg = ncg_mtop(mtop);
2435 snew(ns->simple_aaj,2*ncg);
2436 for(jcg=0; (jcg<ncg); jcg++)
2438 ns->simple_aaj[jcg] = jcg;
2439 ns->simple_aaj[jcg+ncg] = jcg;
2443 /* Create array that determines whether or not atoms have VdW */
2444 snew(ns->bHaveVdW,fr->ntype);
2445 for(i=0; (i<fr->ntype); i++)
2447 for(j=0; (j<fr->ntype); j++)
2449 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2451 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2452 (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2453 (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2454 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2455 (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2459 pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2463 if (!DOMAINDECOMP(cr))
2465 /* This could be reduced with particle decomposition */
2466 ns_realloc_natoms(ns,mtop->natoms);
2469 ns->nblist_initialized=FALSE;
2471 /* nbr list debug dump */
2473 char *ptr=getenv("GMX_DUMP_NL");
2476 ns->dump_nl=strtol(ptr,NULL,10);
2479 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2490 int search_neighbours(FILE *log,t_forcerec *fr,
2491 rvec x[],matrix box,
2492 gmx_localtop_t *top,
2493 gmx_groups_t *groups,
2495 t_nrnb *nrnb,t_mdatoms *md,
2496 real lambda,real *dvdlambda,
2497 gmx_grppairener_t *grppener,
2499 gmx_bool bDoLongRange,
2500 gmx_bool bDoForces,rvec *f)
2502 t_block *cgs=&(top->cgs);
2503 rvec box_size,grid_x0,grid_x1;
2505 real min_size,grid_dens;
2509 gmx_bool *i_egp_flags;
2510 int cg_start,cg_end,start,end;
2513 gmx_domdec_zones_t *dd_zones;
2514 put_in_list_t *put_in_list;
2518 /* Set some local variables */
2520 ngid = groups->grps[egcENER].nr;
2522 for(m=0; (m<DIM); m++)
2524 box_size[m] = box[m][m];
2527 if (fr->ePBC != epbcNONE)
2529 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2531 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.");
2535 min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2536 if (2*fr->rlistlong >= min_size)
2537 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2541 if (DOMAINDECOMP(cr))
2543 ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2547 /* Reset the neighbourlists */
2548 reset_neighbor_list(fr,FALSE,-1,-1);
2550 if (bGrid && bFillGrid)
2554 if (DOMAINDECOMP(cr))
2556 dd_zones = domdec_zones(cr->dd);
2562 get_nsgrid_boundaries(grid,NULL,box,NULL,NULL,NULL,
2563 cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2565 grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2566 fr->rlistlong,grid_dens);
2570 /* Don't know why this all is... (DvdS 3/99) */
2576 end = (cgs->nr+1)/2;
2579 if (DOMAINDECOMP(cr))
2582 fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2584 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2588 fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2589 grid->icg0 = fr->cg0;
2590 grid->icg1 = fr->hcg;
2598 calc_elemnr(log,grid,start,end,cgs->nr);
2600 grid_last(log,grid,start,end,cgs->nr);
2604 check_grid(debug,grid);
2605 print_grid(debug,grid);
2610 /* Set the grid cell index for the test particle only.
2611 * The cell to cg index is not corrected, but that does not matter.
2613 fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2617 if (!fr->ns.bCGlist)
2619 put_in_list = put_in_list_at;
2623 put_in_list = put_in_list_cg;
2630 nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2631 grid,x,ns->bexcl,ns->bExcludeAlleg,
2632 nrnb,md,lambda,dvdlambda,grppener,
2633 put_in_list,ns->bHaveVdW,
2634 bDoLongRange,bDoForces,f,
2637 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2638 * the classical calculation. The charge-charge interaction
2639 * between QM and MM atoms is handled in the QMMM core calculation
2640 * (see QMMM.c). The VDW however, we'd like to compute classically
2641 * and the QM MM atom pairs have just been put in the
2642 * corresponding neighbourlists. in case of QMMM we still need to
2643 * fill a special QMMM neighbourlist that contains all neighbours
2644 * of the QM atoms. If bQMMM is true, this list will now be made:
2646 if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
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_qmmm,ns->bHaveVdW,
2652 bDoLongRange,bDoForces,f,
2658 nsearch = ns_simple_core(fr,top,md,box,box_size,
2659 ns->bexcl,ns->simple_aaj,
2660 ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2668 inc_nrnb(nrnb,eNR_NS,nsearch);
2669 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2674 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2675 matrix scale_tot,rvec *x)
2677 int cg0,cg1,cg,a0,a1,a,i,j;
2678 real rint,hbuf2,scale;
2680 gmx_bool bIsotropic;
2685 rint = max(ir->rcoulomb,ir->rvdw);
2686 if (ir->rlist < rint)
2688 gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2696 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2698 hbuf2 = sqr(0.5*(ir->rlist - rint));
2699 for(cg=cg0; cg<cg1; cg++)
2701 a0 = cgs->index[cg];
2702 a1 = cgs->index[cg+1];
2703 for(a=a0; a<a1; a++)
2705 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2715 scale = scale_tot[0][0];
2716 for(i=1; i<DIM; i++)
2718 /* With anisotropic scaling, the original spherical ns volumes become
2719 * ellipsoids. To avoid costly transformations we use the minimum
2720 * eigenvalue of the scaling matrix for determining the buffer size.
2721 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2723 scale = min(scale,scale_tot[i][i]);
2724 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2730 if (scale_tot[i][j] != 0)
2736 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2739 for(cg=cg0; cg<cg1; cg++)
2741 svmul(scale,cg_cm[cg],cgsc);
2742 a0 = cgs->index[cg];
2743 a1 = cgs->index[cg+1];
2744 for(a=a0; a<a1; a++)
2746 if (distance2(cgsc,x[a]) > hbuf2)
2755 /* Anistropic scaling */
2756 for(cg=cg0; cg<cg1; cg++)
2758 /* Since scale_tot contains the transpose of the scaling matrix,
2759 * we need to multiply with the transpose.
2761 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2762 a0 = cgs->index[cg];
2763 a1 = cgs->index[cg+1];
2764 for(a=a0; a<a1; a++)
2766 if (distance2(cgsc,x[a]) > hbuf2)