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.
481 * Avoid the creation of i-lists with no j-particles.
485 /* Some assembly kernels do not support empty lists,
486 * make sure here that we don't generate any empty lists.
487 * With the current ns code this branch is taken in two cases:
488 * No i-particles at all: nri=-1 here
489 * There are i-particles, but no j-particles; nri=0 here
495 /* Close list number nri by incrementing the count */
500 static inline void close_neighbor_list(t_forcerec *fr,gmx_bool bLR,int nls,int eNL,
501 gmx_bool bMakeQMMMnblist)
505 if (bMakeQMMMnblist) {
508 close_nblist(&(fr->QMMMlist));
515 close_nblist(&(fr->nblists[nls].nlist_lr[eNL]));
519 for(n=0; n<fr->nnblists; n++)
521 for(i=0; (i<eNL_NR); i++)
523 close_nblist(&(fr->nblists[n].nlist_sr[i]));
530 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
534 if (nlist->nrj >= nlist->maxnrj)
536 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
538 fprintf(debug,"Increasing %s nblist %s j size to %d\n",
539 bLR ? "LR" : "SR",nrnb_str(nlist->il_code),nlist->maxnrj);
541 srenew(nlist->jjnr,nlist->maxnrj);
544 nlist->jjnr[nrj] = j_atom;
548 static inline void add_j_to_nblist_cg(t_nblist *nlist,
549 atom_id j_start,int j_end,
550 t_excl *bexcl,gmx_bool bLR)
555 if (nlist->nrj >= nlist->maxnrj)
557 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
559 fprintf(debug,"Increasing %s nblist %s j size to %d\n",
560 bLR ? "LR" : "SR",nrnb_str(nlist->il_code),nlist->maxnrj);
562 srenew(nlist->jjnr ,nlist->maxnrj);
563 srenew(nlist->jjnr_end,nlist->maxnrj);
564 srenew(nlist->excl ,nlist->maxnrj*MAX_CGCGSIZE);
567 nlist->jjnr[nrj] = j_start;
568 nlist->jjnr_end[nrj] = j_end;
570 if (j_end - j_start > MAX_CGCGSIZE)
572 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);
575 /* Set the exclusions */
576 for(j=j_start; j<j_end; j++)
578 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
585 put_in_list_t(gmx_bool bHaveVdW[],
601 put_in_list_at(gmx_bool bHaveVdW[],
616 /* The a[] index has been removed,
617 * to put it back in i_atom should be a[i0] and jj should be a[jj].
622 t_nblist * vdwc_free = NULL;
623 t_nblist * vdw_free = NULL;
624 t_nblist * coul_free = NULL;
625 t_nblist * vdwc_ww = NULL;
626 t_nblist * coul_ww = NULL;
628 int i,j,jcg,igid,gid,nbl_ind,ind_ij;
629 atom_id jj,jj0,jj1,i_atom;
634 real *charge,*chargeB;
636 gmx_bool bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
637 gmx_bool bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
641 /* Copy some pointers */
643 charge = md->chargeA;
644 chargeB = md->chargeB;
647 bPert = md->bPerturbed;
651 nicg = index[icg+1]-i0;
653 /* Get the i charge group info */
654 igid = GET_CGINFO_GID(cginfo[icg]);
655 iwater = GET_CGINFO_SOLOPT(cginfo[icg]);
660 /* Check if any of the particles involved are perturbed.
661 * If not we can do the cheaper normal put_in_list
662 * and use more solvent optimization.
664 for(i=0; i<nicg; i++)
666 bFreeEnergy |= bPert[i0+i];
668 /* Loop over the j charge groups */
669 for(j=0; (j<nj && !bFreeEnergy); j++)
674 /* Finally loop over the atoms in the j-charge group */
675 for(jj=jj0; jj<jj1; jj++)
677 bFreeEnergy |= bPert[jj];
682 /* Unpack pointers to neighbourlist structs */
683 if (fr->nnblists == 1)
689 nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
693 nlist = fr->nblists[nbl_ind].nlist_lr;
697 nlist = fr->nblists[nbl_ind].nlist_sr;
700 if (iwater != esolNO)
702 vdwc = &nlist[eNL_VDWQQ_WATER];
703 vdw = &nlist[eNL_VDW];
704 coul = &nlist[eNL_QQ_WATER];
705 #ifndef DISABLE_WATERWATER_NLIST
706 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
707 coul_ww = &nlist[eNL_QQ_WATERWATER];
712 vdwc = &nlist[eNL_VDWQQ];
713 vdw = &nlist[eNL_VDW];
714 coul = &nlist[eNL_QQ];
719 if (iwater != esolNO)
721 /* Loop over the atoms in the i charge group */
723 gid = GID(igid,jgid,ngid);
724 /* Create new i_atom for each energy group */
725 if (bDoCoul && bDoVdW)
727 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
728 #ifndef DISABLE_WATERWATER_NLIST
729 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
734 new_i_nblist(vdw,bLR,i_atom,shift,gid);
738 new_i_nblist(coul,bLR,i_atom,shift,gid);
739 #ifndef DISABLE_WATERWATER_NLIST
740 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
743 /* Loop over the j charge groups */
744 for(j=0; (j<nj); j++)
754 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
756 if (iwater == esolSPC && jwater == esolSPC)
758 /* Interaction between two SPC molecules */
761 /* VdW only - only first atoms in each water interact */
762 add_j_to_nblist(vdw,jj0,bLR);
766 #ifdef DISABLE_WATERWATER_NLIST
767 /* Add entries for the three atoms - only do VdW if we need to */
770 add_j_to_nblist(coul,jj0,bLR);
774 add_j_to_nblist(vdwc,jj0,bLR);
776 add_j_to_nblist(coul,jj0+1,bLR);
777 add_j_to_nblist(coul,jj0+2,bLR);
779 /* One entry for the entire water-water interaction */
782 add_j_to_nblist(coul_ww,jj0,bLR);
786 add_j_to_nblist(vdwc_ww,jj0,bLR);
791 else if (iwater == esolTIP4P && jwater == esolTIP4P)
793 /* Interaction between two TIP4p molecules */
796 /* VdW only - only first atoms in each water interact */
797 add_j_to_nblist(vdw,jj0,bLR);
801 #ifdef DISABLE_WATERWATER_NLIST
802 /* Add entries for the four atoms - only do VdW if we need to */
805 add_j_to_nblist(vdw,jj0,bLR);
807 add_j_to_nblist(coul,jj0+1,bLR);
808 add_j_to_nblist(coul,jj0+2,bLR);
809 add_j_to_nblist(coul,jj0+3,bLR);
811 /* One entry for the entire water-water interaction */
814 add_j_to_nblist(coul_ww,jj0,bLR);
818 add_j_to_nblist(vdwc_ww,jj0,bLR);
825 /* j charge group is not water, but i is.
826 * Add entries to the water-other_atom lists; the geometry of the water
827 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
828 * so we don't care if it is SPC or TIP4P...
835 for(jj=jj0; (jj<jj1); jj++)
839 add_j_to_nblist(coul,jj,bLR);
845 for(jj=jj0; (jj<jj1); jj++)
847 if (bHaveVdW[type[jj]])
849 add_j_to_nblist(vdw,jj,bLR);
855 /* _charge_ _groups_ interact with both coulomb and LJ */
856 /* Check which atoms we should add to the lists! */
857 for(jj=jj0; (jj<jj1); jj++)
859 if (bHaveVdW[type[jj]])
863 add_j_to_nblist(vdwc,jj,bLR);
867 add_j_to_nblist(vdw,jj,bLR);
870 else if (charge[jj] != 0)
872 add_j_to_nblist(coul,jj,bLR);
879 close_i_nblist(coul);
880 close_i_nblist(vdwc);
881 #ifndef DISABLE_WATERWATER_NLIST
882 close_i_nblist(coul_ww);
883 close_i_nblist(vdwc_ww);
888 /* no solvent as i charge group */
889 /* Loop over the atoms in the i charge group */
890 for(i=0; i<nicg; i++)
893 gid = GID(igid,jgid,ngid);
896 /* Create new i_atom for each energy group */
897 if (bDoVdW && bDoCoul)
899 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
903 new_i_nblist(vdw,bLR,i_atom,shift,gid);
907 new_i_nblist(coul,bLR,i_atom,shift,gid);
909 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
910 bDoCoul_i = (bDoCoul && qi!=0);
912 if (bDoVdW_i || bDoCoul_i)
914 /* Loop over the j charge groups */
915 for(j=0; (j<nj); j++)
919 /* Check for large charge groups */
930 /* Finally loop over the atoms in the j-charge group */
931 for(jj=jj0; jj<jj1; jj++)
933 bNotEx = NOTEXCL(bExcl,i,jj);
941 add_j_to_nblist(coul,jj,bLR);
946 if (bHaveVdW[type[jj]])
948 add_j_to_nblist(vdw,jj,bLR);
953 if (bHaveVdW[type[jj]])
957 add_j_to_nblist(vdwc,jj,bLR);
961 add_j_to_nblist(vdw,jj,bLR);
964 else if (charge[jj] != 0)
966 add_j_to_nblist(coul,jj,bLR);
974 close_i_nblist(coul);
975 close_i_nblist(vdwc);
981 /* we are doing free energy */
982 vdwc_free = &nlist[eNL_VDWQQ_FREE];
983 vdw_free = &nlist[eNL_VDW_FREE];
984 coul_free = &nlist[eNL_QQ_FREE];
985 /* Loop over the atoms in the i charge group */
986 for(i=0; i<nicg; i++)
989 gid = GID(igid,jgid,ngid);
991 qiB = chargeB[i_atom];
993 /* Create new i_atom for each energy group */
994 if (bDoVdW && bDoCoul)
995 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
997 new_i_nblist(vdw,bLR,i_atom,shift,gid);
999 new_i_nblist(coul,bLR,i_atom,shift,gid);
1001 new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
1002 new_i_nblist(coul_free,bLR,i_atom,shift,gid);
1003 new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
1005 bDoVdW_i = (bDoVdW &&
1006 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
1007 bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
1008 /* For TIP4P the first atom does not have a charge,
1009 * but the last three do. So we should still put an atom
1010 * without LJ but with charge in the water-atom neighborlist
1011 * for a TIP4p i charge group.
1012 * For SPC type water the first atom has LJ and charge,
1013 * so there is no such problem.
1015 if (iwater == esolNO)
1017 bDoCoul_i_sol = bDoCoul_i;
1021 bDoCoul_i_sol = bDoCoul;
1024 if (bDoVdW_i || bDoCoul_i_sol)
1026 /* Loop over the j charge groups */
1027 for(j=0; (j<nj); j++)
1031 /* Check for large charge groups */
1042 /* Finally loop over the atoms in the j-charge group */
1043 bFree = bPert[i_atom];
1044 for(jj=jj0; (jj<jj1); jj++)
1046 bFreeJ = bFree || bPert[jj];
1047 /* Complicated if, because the water H's should also
1048 * see perturbed j-particles
1050 if (iwater==esolNO || i==0 || bFreeJ)
1052 bNotEx = NOTEXCL(bExcl,i,jj);
1060 if (charge[jj]!=0 || chargeB[jj]!=0)
1062 add_j_to_nblist(coul_free,jj,bLR);
1065 else if (!bDoCoul_i)
1067 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1069 add_j_to_nblist(vdw_free,jj,bLR);
1074 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1076 if (charge[jj]!=0 || chargeB[jj]!=0)
1078 add_j_to_nblist(vdwc_free,jj,bLR);
1082 add_j_to_nblist(vdw_free,jj,bLR);
1085 else if (charge[jj]!=0 || chargeB[jj]!=0)
1086 add_j_to_nblist(coul_free,jj,bLR);
1091 /* This is done whether or not bWater is set */
1092 if (charge[jj] != 0)
1094 add_j_to_nblist(coul,jj,bLR);
1097 else if (!bDoCoul_i_sol)
1099 if (bHaveVdW[type[jj]])
1101 add_j_to_nblist(vdw,jj,bLR);
1106 if (bHaveVdW[type[jj]])
1108 if (charge[jj] != 0)
1110 add_j_to_nblist(vdwc,jj,bLR);
1114 add_j_to_nblist(vdw,jj,bLR);
1117 else if (charge[jj] != 0)
1119 add_j_to_nblist(coul,jj,bLR);
1127 close_i_nblist(vdw);
1128 close_i_nblist(coul);
1129 close_i_nblist(vdwc);
1130 close_i_nblist(vdw_free);
1131 close_i_nblist(coul_free);
1132 close_i_nblist(vdwc_free);
1138 put_in_list_qmmm(gmx_bool bHaveVdW[],
1154 int i,j,jcg,igid,gid;
1155 atom_id jj,jj0,jj1,i_atom;
1159 /* Get atom range */
1161 nicg = index[icg+1]-i0;
1163 /* Get the i charge group info */
1164 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1166 coul = &fr->QMMMlist;
1168 /* Loop over atoms in the ith charge group */
1169 for (i=0;i<nicg;i++)
1172 gid = GID(igid,jgid,ngid);
1173 /* Create new i_atom for each energy group */
1174 new_i_nblist(coul,bLR,i_atom,shift,gid);
1176 /* Loop over the j charge groups */
1181 /* Charge groups cannot have QM and MM atoms simultaneously */
1186 /* Finally loop over the atoms in the j-charge group */
1187 for(jj=jj0; jj<jj1; jj++)
1189 bNotEx = NOTEXCL(bExcl,i,jj);
1191 add_j_to_nblist(coul,jj,bLR);
1195 close_i_nblist(coul);
1200 put_in_list_cg(gmx_bool bHaveVdW[],
1216 int igid,gid,nbl_ind;
1220 cginfo = fr->cginfo[icg];
1222 igid = GET_CGINFO_GID(cginfo);
1223 gid = GID(igid,jgid,ngid);
1225 /* Unpack pointers to neighbourlist structs */
1226 if (fr->nnblists == 1)
1232 nbl_ind = fr->gid2nblists[gid];
1236 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1240 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1243 /* Make a new neighbor list for charge group icg.
1244 * Currently simply one neighbor list is made with LJ and Coulomb.
1245 * If required, zero interactions could be removed here
1246 * or in the force loop.
1248 new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1249 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1251 for(j=0; (j<nj); j++)
1254 /* Skip the icg-icg pairs if all self interactions are excluded */
1255 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1257 /* Here we add the j charge group jcg to the list,
1258 * exclusions are also added to the list.
1260 add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,bLR);
1264 close_i_nblist(vdwc);
1267 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1274 for(i=start; i<end; i++)
1276 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1278 SETEXCL(bexcl,i-start,excl->a[k]);
1284 for(i=start; i<end; i++)
1286 for(k=excl->index[i]; k<excl->index[i+1]; k++)
1288 RMEXCL(bexcl,i-start,excl->a[k]);
1294 int calc_naaj(int icg,int cgtot)
1298 if ((cgtot % 2) == 1)
1300 /* Odd number of charge groups, easy */
1301 naaj = 1 + (cgtot/2);
1303 else if ((cgtot % 4) == 0)
1305 /* Multiple of four is hard */
1342 fprintf(log,"naaj=%d\n",naaj);
1348 /************************************************
1350 * S I M P L E C O R E S T U F F
1352 ************************************************/
1354 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1355 rvec b_inv,int *shift)
1357 /* This code assumes that the cut-off is smaller than
1358 * a half times the smallest diagonal element of the box.
1365 /* Compute diff vector */
1366 dz = xj[ZZ] - xi[ZZ];
1367 dy = xj[YY] - xi[YY];
1368 dx = xj[XX] - xi[XX];
1370 /* Perform NINT operation, using trunc operation, therefore
1371 * we first add 2.5 then subtract 2 again
1373 tz = dz*b_inv[ZZ] + h25;
1375 dz -= tz*box[ZZ][ZZ];
1376 dy -= tz*box[ZZ][YY];
1377 dx -= tz*box[ZZ][XX];
1379 ty = dy*b_inv[YY] + h25;
1381 dy -= ty*box[YY][YY];
1382 dx -= ty*box[YY][XX];
1384 tx = dx*b_inv[XX]+h25;
1386 dx -= tx*box[XX][XX];
1388 /* Distance squared */
1389 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1391 *shift = XYZ2IS(tx,ty,tz);
1396 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1397 rvec b_inv,int *shift)
1405 /* Compute diff vector */
1406 dx = xj[XX] - xi[XX];
1407 dy = xj[YY] - xi[YY];
1408 dz = xj[ZZ] - xi[ZZ];
1410 /* Perform NINT operation, using trunc operation, therefore
1411 * we first add 1.5 then subtract 1 again
1413 tx = dx*b_inv[XX] + h15;
1414 ty = dy*b_inv[YY] + h15;
1415 tz = dz*b_inv[ZZ] + h15;
1420 /* Correct diff vector for translation */
1421 ddx = tx*box_size[XX] - dx;
1422 ddy = ty*box_size[YY] - dy;
1423 ddz = tz*box_size[ZZ] - dz;
1425 /* Distance squared */
1426 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1428 *shift = XYZ2IS(tx,ty,tz);
1433 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1434 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1435 int icg,int jgid,t_block *cgs,t_excl bexcl[],
1436 int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1438 if (nsbuf->nj + nrj > MAX_CG)
1440 put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1441 cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE);
1442 /* Reset buffer contents */
1443 nsbuf->ncg = nsbuf->nj = 0;
1445 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1449 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1450 int njcg,atom_id jcg[],
1451 matrix box,rvec b_inv,real rcut2,
1452 t_block *cgs,t_ns_buf **ns_buf,
1453 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1454 t_excl bexcl[],t_forcerec *fr,
1455 put_in_list_t *put_in_list)
1459 int *cginfo=fr->cginfo;
1460 atom_id cg_j,*cgindex;
1463 cgindex = cgs->index;
1465 for(j=0; (j<njcg); j++)
1468 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1469 if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1471 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1472 if (!(i_egp_flags[jgid] & EGP_EXCL))
1474 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1475 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1482 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1483 int njcg,atom_id jcg[],
1484 gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1485 t_block *cgs,t_ns_buf **ns_buf,
1486 gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1487 t_excl bexcl[],t_forcerec *fr,
1488 put_in_list_t *put_in_list)
1492 int *cginfo=fr->cginfo;
1493 atom_id cg_j,*cgindex;
1496 cgindex = cgs->index;
1500 for(j=0; (j<njcg); j++)
1503 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1504 if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1506 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1507 if (!(i_egp_flags[jgid] & EGP_EXCL))
1509 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1510 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1518 for(j=0; (j<njcg); j++)
1521 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1522 if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1523 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1524 if (!(i_egp_flags[jgid] & EGP_EXCL))
1526 add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1527 bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1535 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1537 static int ns_simple_core(t_forcerec *fr,
1538 gmx_localtop_t *top,
1540 matrix box,rvec box_size,
1541 t_excl bexcl[],atom_id *aaj,
1542 int ngid,t_ns_buf **ns_buf,
1543 put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1547 int nsearch,icg,jcg,igid,i0,nri,nn;
1550 /* atom_id *i_atoms; */
1551 t_block *cgs=&(top->cgs);
1552 t_blocka *excl=&(top->excls);
1555 gmx_bool bBox,bTriclinic;
1558 rlist2 = sqr(fr->rlist);
1560 bBox = (fr->ePBC != epbcNONE);
1563 for(m=0; (m<DIM); m++)
1565 b_inv[m] = divide_err(1.0,box_size[m]);
1567 bTriclinic = TRICLINIC(box);
1574 cginfo = fr->cginfo;
1577 for (icg=fr->cg0; (icg<fr->hcg); icg++)
1580 i0 = cgs->index[icg];
1581 nri = cgs->index[icg+1]-i0;
1582 i_atoms = &(cgs->a[i0]);
1583 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1584 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1586 igid = GET_CGINFO_GID(cginfo[icg]);
1587 i_egp_flags = fr->egp_flags + ngid*igid;
1588 setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1590 naaj=calc_naaj(icg,cgs->nr);
1593 ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1594 box,b_inv,rlist2,cgs,ns_buf,
1595 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1599 ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1600 bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1601 bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1605 for(nn=0; (nn<ngid); nn++)
1607 for(k=0; (k<SHIFTS); k++)
1609 nsbuf = &(ns_buf[nn][k]);
1612 put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1613 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE);
1614 nsbuf->ncg=nsbuf->nj=0;
1618 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1619 setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1621 close_neighbor_list(fr,FALSE,-1,-1,FALSE);
1626 /************************************************
1628 * N S 5 G R I D S T U F F
1630 ************************************************/
1632 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1633 int *dx0,int *dx1,real *dcx2)
1661 for(i=xgi0; i>=0; i--)
1663 dcx = (i+1)*gridx-x;
1670 for(i=xgi1; i<Nx; i++)
1683 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1684 int ncpddc,int shift_min,int shift_max,
1685 int *g0,int *g1,real *dcx2)
1688 int g_min,g_max,shift_home;
1721 g_min = (shift_min == shift_home ? 0 : ncpddc);
1722 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1729 else if (shift_max < 0)
1744 /* Check one grid cell down */
1745 dcx = ((*g0 - 1) + 1)*gridx - x;
1757 /* Check one grid cell up */
1758 dcx = (*g1 + 1)*gridx - x;
1770 #define sqr(x) ((x)*(x))
1771 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1772 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1773 /****************************************************
1775 * F A S T N E I G H B O R S E A R C H I N G
1777 * Optimized neighboursearching routine using grid
1778 * at least 1x1x1, see GROMACS manual
1780 ****************************************************/
1782 static void do_longrange(t_commrec *cr,gmx_localtop_t *top,t_forcerec *fr,
1783 int ngid,t_mdatoms *md,int icg,
1785 atom_id lr[],t_excl bexcl[],int shift,
1786 rvec x[],rvec box_size,t_nrnb *nrnb,
1787 real lambda,real *dvdlambda,
1788 gmx_grppairener_t *grppener,
1789 gmx_bool bDoVdW,gmx_bool bDoCoul,
1790 gmx_bool bEvaluateNow,put_in_list_t *put_in_list,
1791 gmx_bool bHaveVdW[],
1792 gmx_bool bDoForces,rvec *f)
1797 for(n=0; n<fr->nnblists; n++)
1799 for(i=0; (i<eNL_NR); i++)
1801 nl = &fr->nblists[n].nlist_lr[i];
1802 if ((nl->nri > nl->maxnri-32) || bEvaluateNow)
1804 close_neighbor_list(fr,TRUE,n,i,FALSE);
1805 /* Evaluate the energies and forces */
1806 do_nonbonded(cr,fr,x,f,md,NULL,
1807 grppener->ener[fr->bBHAM ? egBHAMLR : egLJLR],
1808 grppener->ener[egCOULLR],
1809 grppener->ener[egGB],box_size,
1810 nrnb,lambda,dvdlambda,n,i,
1811 GMX_DONB_LR | GMX_DONB_FORCES);
1813 reset_neighbor_list(fr,TRUE,n,i);
1820 /* Put the long range particles in a list */
1821 /* do_longrange is never called for QMMM */
1822 put_in_list(bHaveVdW,ngid,md,icg,jgid,nlr,lr,top->cgs.index,
1823 bexcl,shift,fr,TRUE,bDoVdW,bDoCoul);
1827 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1828 real *rvdw2,real *rcoul2,
1829 real *rs2,real *rm2,real *rl2)
1831 *rs2 = sqr(fr->rlist);
1832 if (bDoLongRange && fr->bTwinRange)
1834 /* The VdW and elec. LR cut-off's could be different,
1835 * so we can not simply set them to rlistlong.
1837 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1838 fr->rvdw > fr->rlist)
1840 *rvdw2 = sqr(fr->rlistlong);
1844 *rvdw2 = sqr(fr->rvdw);
1846 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1847 fr->rcoulomb > fr->rlist)
1849 *rcoul2 = sqr(fr->rlistlong);
1853 *rcoul2 = sqr(fr->rcoulomb);
1858 /* Workaround for a gcc -O3 or -ffast-math problem */
1862 *rm2 = min(*rvdw2,*rcoul2);
1863 *rl2 = max(*rvdw2,*rcoul2);
1866 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
1868 real rvdw2,rcoul2,rs2,rm2,rl2;
1871 get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1873 /* Short range buffers */
1874 snew(ns->nl_sr,ngid);
1877 snew(ns->nlr_ljc,ngid);
1878 snew(ns->nlr_one,ngid);
1882 /* Long range VdW and Coul buffers */
1883 snew(ns->nl_lr_ljc,ngid);
1887 /* Long range VdW or Coul only buffers */
1888 snew(ns->nl_lr_one,ngid);
1890 for(j=0; (j<ngid); j++) {
1891 snew(ns->nl_sr[j],MAX_CG);
1894 snew(ns->nl_lr_ljc[j],MAX_CG);
1898 snew(ns->nl_lr_one[j],MAX_CG);
1904 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
1909 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
1910 matrix box,rvec box_size,int ngid,
1911 gmx_localtop_t *top,
1912 t_grid *grid,rvec x[],
1913 t_excl bexcl[],gmx_bool *bExcludeAlleg,
1914 t_nrnb *nrnb,t_mdatoms *md,
1915 real lambda,real *dvdlambda,
1916 gmx_grppairener_t *grppener,
1917 put_in_list_t *put_in_list,
1918 gmx_bool bHaveVdW[],
1919 gmx_bool bDoLongRange,gmx_bool bDoForces,rvec *f,
1920 gmx_bool bMakeQMMMnblist)
1923 atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
1924 int *nlr_ljc,*nlr_one,*nsr;
1925 gmx_domdec_t *dd=NULL;
1926 t_block *cgs=&(top->cgs);
1927 int *cginfo=fr->cginfo;
1928 /* atom_id *i_atoms,*cgsindex=cgs->index; */
1930 int cell_x,cell_y,cell_z;
1931 int d,tx,ty,tz,dx,dy,dz,cj;
1932 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
1933 int zsh_ty,zsh_tx,ysh_tx;
1935 int dx0,dx1,dy0,dy1,dz0,dz1;
1936 int Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
1937 real gridx,gridy,gridz,grid_x,grid_y,grid_z;
1938 real *dcx2,*dcy2,*dcz2;
1940 int cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
1941 int jcg0,jcg1,jjcg,cgj0,jgid;
1942 int *grida,*gridnra,*gridind;
1943 gmx_bool rvdw_lt_rcoul,rcoul_lt_rvdw;
1944 rvec xi,*cgcm,grid_offset;
1945 real r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
1947 gmx_bool bDomDec,bTriclinicX,bTriclinicY;
1952 bDomDec = DOMAINDECOMP(cr);
1958 bTriclinicX = ((YY < grid->npbcdim &&
1959 (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
1960 (ZZ < grid->npbcdim &&
1961 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
1962 bTriclinicY = (ZZ < grid->npbcdim &&
1963 (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
1967 get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
1969 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
1970 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
1972 if (bMakeQMMMnblist)
1980 nl_lr_ljc = ns->nl_lr_ljc;
1981 nl_lr_one = ns->nl_lr_one;
1982 nlr_ljc = ns->nlr_ljc;
1983 nlr_one = ns->nlr_one;
1991 gridind = grid->index;
1992 gridnra = grid->nra;
1995 gridx = grid->cell_size[XX];
1996 gridy = grid->cell_size[YY];
1997 gridz = grid->cell_size[ZZ];
2001 copy_rvec(grid->cell_offset,grid_offset);
2002 copy_ivec(grid->ncpddc,ncpddc);
2007 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2008 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2009 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2010 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2011 if (zsh_tx!=0 && ysh_tx!=0)
2013 /* This could happen due to rounding, when both ratios are 0.5 */
2022 /* We only want a list for the test particle */
2031 /* Set the shift range */
2032 for(d=0; d<DIM; d++)
2036 /* Check if we need periodicity shifts.
2037 * Without PBC or with domain decomposition we don't need them.
2039 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2046 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2057 /* Loop over charge groups */
2058 for(icg=cg0; (icg < cg1); icg++)
2060 igid = GET_CGINFO_GID(cginfo[icg]);
2061 /* Skip this charge group if all energy groups are excluded! */
2062 if (bExcludeAlleg[igid])
2067 i0 = cgs->index[icg];
2069 if (bMakeQMMMnblist)
2071 /* Skip this charge group if it is not a QM atom while making a
2072 * QM/MM neighbourlist
2074 if (md->bQM[i0]==FALSE)
2076 continue; /* MM particle, go to next particle */
2079 /* Compute the number of charge groups that fall within the control
2082 naaj = calc_naaj(icg,cgsnr);
2089 /* make a normal neighbourlist */
2093 /* Get the j charge-group and dd cell shift ranges */
2094 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
2099 /* Compute the number of charge groups that fall within the control
2102 naaj = calc_naaj(icg,cgsnr);
2108 /* The i-particle is awlways the test particle,
2109 * so we want all j-particles
2111 max_jcg = cgsnr - 1;
2115 max_jcg = jcg1 - cgsnr;
2120 i_egp_flags = fr->egp_flags + igid*ngid;
2122 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2123 setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
2125 ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2127 /* Changed iicg to icg, DvdS 990115
2128 * (but see consistency check above, DvdS 990330)
2131 fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2132 icg,naaj,cell_x,cell_y,cell_z);
2134 /* Loop over shift vectors in three dimensions */
2135 for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2137 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2138 /* Calculate range of cells in Z direction that have the shift tz */
2139 zgi = cell_z + tz*Nz;
2142 get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2144 get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2145 ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2151 for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2153 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2154 /* Calculate range of cells in Y direction that have the shift ty */
2157 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2161 ygi = cell_y + ty*Ny;
2164 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2166 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2167 ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2173 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2175 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2176 /* Calculate range of cells in X direction that have the shift tx */
2179 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2183 xgi = cell_x + tx*Nx;
2186 get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2188 get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2189 ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2195 /* Get shift vector */
2196 shift=XYZ2IS(tx,ty,tz);
2198 range_check(shift,0,SHIFTS);
2200 for(nn=0; (nn<ngid); nn++)
2207 fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2208 shift,dx0,dx1,dy0,dy1,dz0,dz1);
2209 fprintf(log,"cgcm: %8.3f %8.3f %8.3f\n",cgcm[icg][XX],
2210 cgcm[icg][YY],cgcm[icg][ZZ]);
2211 fprintf(log,"xi: %8.3f %8.3f %8.3f\n",XI,YI,ZI);
2213 for (dx=dx0; (dx<=dx1); dx++)
2215 tmp1 = rl2 - dcx2[dx];
2216 for (dy=dy0; (dy<=dy1); dy++)
2218 tmp2 = tmp1 - dcy2[dy];
2221 for (dz=dz0; (dz<=dz1); dz++) {
2222 if (tmp2 > dcz2[dz]) {
2223 /* Find grid-cell cj in which possible neighbours are */
2224 cj = xyz2ci(Ny,Nz,dx,dy,dz);
2226 /* Check out how many cgs (nrj) there in this cell */
2229 /* Find the offset in the cg list */
2232 /* Check if all j's are out of range so we
2233 * can skip the whole cell.
2234 * Should save some time, especially with DD.
2237 (grida[cgj0] >= max_jcg &&
2238 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2244 for (j=0; (j<nrj); j++)
2246 jjcg = grida[cgj0+j];
2248 /* check whether this guy is in range! */
2249 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2252 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2254 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2255 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2256 /* check energy group exclusions */
2257 if (!(i_egp_flags[jgid] & EGP_EXCL))
2261 if (nsr[jgid] >= MAX_CG)
2263 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2264 nsr[jgid],nl_sr[jgid],
2265 cgs->index,/* cgsatoms, */ bexcl,
2266 shift,fr,FALSE,TRUE,TRUE);
2269 nl_sr[jgid][nsr[jgid]++]=jjcg;
2273 if (nlr_ljc[jgid] >= MAX_CG)
2275 do_longrange(cr,top,fr,ngid,md,icg,jgid,
2277 nl_lr_ljc[jgid],bexcl,shift,x,
2287 nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2291 if (nlr_one[jgid] >= MAX_CG) {
2292 do_longrange(cr,top,fr,ngid,md,icg,jgid,
2294 nl_lr_one[jgid],bexcl,shift,x,
2298 rvdw_lt_rcoul,rcoul_lt_rvdw,FALSE,
2304 nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2316 /* CHECK whether there is anything left in the buffers */
2317 for(nn=0; (nn<ngid); nn++)
2321 put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2322 cgs->index, /* cgsatoms, */ bexcl,
2323 shift,fr,FALSE,TRUE,TRUE);
2326 if (nlr_ljc[nn] > 0)
2328 do_longrange(cr,top,fr,ngid,md,icg,nn,nlr_ljc[nn],
2329 nl_lr_ljc[nn],bexcl,shift,x,box_size,nrnb,
2330 lambda,dvdlambda,grppener,TRUE,TRUE,FALSE,
2331 put_in_list,bHaveVdW,bDoForces,f);
2334 if (nlr_one[nn] > 0)
2336 do_longrange(cr,top,fr,ngid,md,icg,nn,nlr_one[nn],
2337 nl_lr_one[nn],bexcl,shift,x,box_size,nrnb,
2338 lambda,dvdlambda,grppener,
2339 rvdw_lt_rcoul,rcoul_lt_rvdw,FALSE,
2340 put_in_list,bHaveVdW,bDoForces,f);
2346 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2347 setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2349 /* Perform any left over force calculations */
2350 for (nn=0; (nn<ngid); nn++)
2354 do_longrange(cr,top,fr,0,md,icg,nn,nlr_ljc[nn],
2355 nl_lr_ljc[nn],bexcl,shift,x,box_size,nrnb,
2356 lambda,dvdlambda,grppener,
2357 TRUE,TRUE,TRUE,put_in_list,bHaveVdW,bDoForces,f);
2360 do_longrange(cr,top,fr,0,md,icg,nn,nlr_one[nn],
2361 nl_lr_one[nn],bexcl,shift,x,box_size,nrnb,
2362 lambda,dvdlambda,grppener,
2363 rvdw_lt_rcoul,rcoul_lt_rvdw,
2364 TRUE,put_in_list,bHaveVdW,bDoForces,f);
2369 /* Close off short range neighbourlists */
2370 close_neighbor_list(fr,FALSE,-1,-1,bMakeQMMMnblist);
2375 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2379 if (natoms > ns->nra_alloc)
2381 ns->nra_alloc = over_alloc_dd(natoms);
2382 srenew(ns->bexcl,ns->nra_alloc);
2383 for(i=0; i<ns->nra_alloc; i++)
2390 void init_ns(FILE *fplog,const t_commrec *cr,
2391 gmx_ns_t *ns,t_forcerec *fr,
2392 const gmx_mtop_t *mtop,
2395 int mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2399 /* Compute largest charge groups size (# atoms) */
2401 for(mt=0; mt<mtop->nmoltype; mt++) {
2402 cgs = &mtop->moltype[mt].cgs;
2403 for (icg=0; (icg < cgs->nr); icg++)
2405 nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2409 /* Verify whether largest charge group is <= max cg.
2410 * This is determined by the type of the local exclusion type
2411 * Exclusions are stored in bits. (If the type is not large
2412 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2414 maxcg = sizeof(t_excl)*8;
2415 if (nr_in_cg > maxcg)
2417 gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2421 ngid = mtop->groups.grps[egcENER].nr;
2422 snew(ns->bExcludeAlleg,ngid);
2423 for(i=0; i<ngid; i++) {
2424 ns->bExcludeAlleg[i] = TRUE;
2425 for(j=0; j<ngid; j++)
2427 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2429 ns->bExcludeAlleg[i] = FALSE;
2436 ns->grid = init_grid(fplog,fr);
2437 init_nsgrid_lists(fr,ngid,ns);
2442 snew(ns->ns_buf,ngid);
2443 for(i=0; (i<ngid); i++)
2445 snew(ns->ns_buf[i],SHIFTS);
2447 ncg = ncg_mtop(mtop);
2448 snew(ns->simple_aaj,2*ncg);
2449 for(jcg=0; (jcg<ncg); jcg++)
2451 ns->simple_aaj[jcg] = jcg;
2452 ns->simple_aaj[jcg+ncg] = jcg;
2456 /* Create array that determines whether or not atoms have VdW */
2457 snew(ns->bHaveVdW,fr->ntype);
2458 for(i=0; (i<fr->ntype); i++)
2460 for(j=0; (j<fr->ntype); j++)
2462 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2464 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2465 (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2466 (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2467 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2468 (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2472 pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2476 if (!DOMAINDECOMP(cr))
2478 /* This could be reduced with particle decomposition */
2479 ns_realloc_natoms(ns,mtop->natoms);
2482 ns->nblist_initialized=FALSE;
2484 /* nbr list debug dump */
2486 char *ptr=getenv("GMX_DUMP_NL");
2489 ns->dump_nl=strtol(ptr,NULL,10);
2492 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2503 int search_neighbours(FILE *log,t_forcerec *fr,
2504 rvec x[],matrix box,
2505 gmx_localtop_t *top,
2506 gmx_groups_t *groups,
2508 t_nrnb *nrnb,t_mdatoms *md,
2509 real lambda,real *dvdlambda,
2510 gmx_grppairener_t *grppener,
2512 gmx_bool bDoLongRange,
2513 gmx_bool bDoForces,rvec *f)
2515 t_block *cgs=&(top->cgs);
2516 rvec box_size,grid_x0,grid_x1;
2518 real min_size,grid_dens;
2522 gmx_bool *i_egp_flags;
2523 int cg_start,cg_end,start,end;
2526 gmx_domdec_zones_t *dd_zones;
2527 put_in_list_t *put_in_list;
2531 /* Set some local variables */
2533 ngid = groups->grps[egcENER].nr;
2535 for(m=0; (m<DIM); m++)
2537 box_size[m] = box[m][m];
2540 if (fr->ePBC != epbcNONE)
2542 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2544 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.");
2548 min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2549 if (2*fr->rlistlong >= min_size)
2550 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2554 if (DOMAINDECOMP(cr))
2556 ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2560 /* Reset the neighbourlists */
2561 reset_neighbor_list(fr,FALSE,-1,-1);
2563 if (bGrid && bFillGrid)
2567 if (DOMAINDECOMP(cr))
2569 dd_zones = domdec_zones(cr->dd);
2575 get_nsgrid_boundaries(grid,NULL,box,NULL,NULL,NULL,
2576 cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2578 grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2579 fr->rlistlong,grid_dens);
2583 /* Don't know why this all is... (DvdS 3/99) */
2589 end = (cgs->nr+1)/2;
2592 if (DOMAINDECOMP(cr))
2595 fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2597 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2601 fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2602 grid->icg0 = fr->cg0;
2603 grid->icg1 = fr->hcg;
2611 calc_elemnr(log,grid,start,end,cgs->nr);
2613 grid_last(log,grid,start,end,cgs->nr);
2617 check_grid(debug,grid);
2618 print_grid(debug,grid);
2623 /* Set the grid cell index for the test particle only.
2624 * The cell to cg index is not corrected, but that does not matter.
2626 fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2630 if (!fr->ns.bCGlist)
2632 put_in_list = put_in_list_at;
2636 put_in_list = put_in_list_cg;
2643 nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2644 grid,x,ns->bexcl,ns->bExcludeAlleg,
2645 nrnb,md,lambda,dvdlambda,grppener,
2646 put_in_list,ns->bHaveVdW,
2647 bDoLongRange,bDoForces,f,
2650 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2651 * the classical calculation. The charge-charge interaction
2652 * between QM and MM atoms is handled in the QMMM core calculation
2653 * (see QMMM.c). The VDW however, we'd like to compute classically
2654 * and the QM MM atom pairs have just been put in the
2655 * corresponding neighbourlists. in case of QMMM we still need to
2656 * fill a special QMMM neighbourlist that contains all neighbours
2657 * of the QM atoms. If bQMMM is true, this list will now be made:
2659 if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
2661 nsearch += nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2662 grid,x,ns->bexcl,ns->bExcludeAlleg,
2663 nrnb,md,lambda,dvdlambda,grppener,
2664 put_in_list_qmmm,ns->bHaveVdW,
2665 bDoLongRange,bDoForces,f,
2671 nsearch = ns_simple_core(fr,top,md,box,box_size,
2672 ns->bexcl,ns->simple_aaj,
2673 ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2681 inc_nrnb(nrnb,eNR_NS,nsearch);
2682 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2687 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2688 matrix scale_tot,rvec *x)
2690 int cg0,cg1,cg,a0,a1,a,i,j;
2691 real rint,hbuf2,scale;
2693 gmx_bool bIsotropic;
2698 rint = max(ir->rcoulomb,ir->rvdw);
2699 if (ir->rlist < rint)
2701 gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2709 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2711 hbuf2 = sqr(0.5*(ir->rlist - rint));
2712 for(cg=cg0; cg<cg1; cg++)
2714 a0 = cgs->index[cg];
2715 a1 = cgs->index[cg+1];
2716 for(a=a0; a<a1; a++)
2718 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2728 scale = scale_tot[0][0];
2729 for(i=1; i<DIM; i++)
2731 /* With anisotropic scaling, the original spherical ns volumes become
2732 * ellipsoids. To avoid costly transformations we use the minimum
2733 * eigenvalue of the scaling matrix for determining the buffer size.
2734 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2736 scale = min(scale,scale_tot[i][i]);
2737 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2743 if (scale_tot[i][j] != 0)
2749 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2752 for(cg=cg0; cg<cg1; cg++)
2754 svmul(scale,cg_cm[cg],cgsc);
2755 a0 = cgs->index[cg];
2756 a1 = cgs->index[cg+1];
2757 for(a=a0; a<a1; a++)
2759 if (distance2(cgsc,x[a]) > hbuf2)
2768 /* Anistropic scaling */
2769 for(cg=cg0; cg<cg1; cg++)
2771 /* Since scale_tot contains the transpose of the scaling matrix,
2772 * we need to multiply with the transpose.
2774 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2775 a0 = cgs->index[cg];
2776 a1 = cgs->index[cg+1];
2777 for(a=a0; a<a1; a++)
2779 if (distance2(cgsc,x[a]) > hbuf2)