1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
50 #include "nonbonded.h"
54 #include "gmx_fatal.h"
57 #include "mtop_util.h"
64 * E X C L U S I O N H A N D L I N G
68 static void SETEXCL_(t_excl e[], atom_id i, atom_id j)
72 static void RMEXCL_(t_excl e[], atom_id i, atom_id j)
74 e[j] = e[j] & ~(1<<i);
76 static gmx_bool ISEXCL_(t_excl e[], atom_id i, atom_id j)
78 return (gmx_bool)(e[j] & (1<<i));
80 static gmx_bool NOTEXCL_(t_excl e[], atom_id i, atom_id j)
82 return !(ISEXCL(e, i, j));
85 #define SETEXCL(e, i, j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
86 #define RMEXCL(e, i, j) (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
87 #define ISEXCL(e, i, j) (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
88 #define NOTEXCL(e, i, j) !(ISEXCL(e, i, j))
92 round_up_to_simd_width(int length, int simd_width)
94 int offset, newlength;
96 offset = (simd_width > 0) ? length % simd_width : 0;
98 return (offset == 0) ? length : length-offset+simd_width;
100 /************************************************
102 * U T I L I T I E S F O R N S
104 ************************************************/
106 static void reallocate_nblist(t_nblist *nl)
110 fprintf(debug, "reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
111 nl->ielec, nl->ivdw, nl->igeometry, nl->type, nl->maxnri);
113 srenew(nl->iinr, nl->maxnri);
114 if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
116 srenew(nl->iinr_end, nl->maxnri);
118 srenew(nl->gid, nl->maxnri);
119 srenew(nl->shift, nl->maxnri);
120 srenew(nl->jindex, nl->maxnri+1);
124 static void init_nblist(FILE *log, t_nblist *nl_sr, t_nblist *nl_lr,
125 int maxsr, int maxlr,
126 int ivdw, int ivdwmod,
127 int ielec, int ielecmod,
128 int igeometry, int type)
134 for (i = 0; (i < 2); i++)
136 nl = (i == 0) ? nl_sr : nl_lr;
137 homenr = (i == 0) ? maxsr : maxlr;
145 /* Set coul/vdw in neighborlist, and for the normal loops we determine
146 * an index of which one to call.
149 nl->ivdwmod = ivdwmod;
151 nl->ielecmod = ielecmod;
153 nl->igeometry = igeometry;
155 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
157 nl->igeometry = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
160 /* This will also set the simd_padding_width field */
161 gmx_nonbonded_set_kernel_pointers( (i == 0) ? log : NULL, nl);
163 /* maxnri is influenced by the number of shifts (maximum is 8)
164 * and the number of energy groups.
165 * If it is not enough, nl memory will be reallocated during the run.
166 * 4 seems to be a reasonable factor, which only causes reallocation
167 * during runs with tiny and many energygroups.
169 nl->maxnri = homenr*4;
178 reallocate_nblist(nl);
183 fprintf(debug, "Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
184 nl->ielec, nl->ivdw, nl->type, gmx_nblist_geometry_names[nl->igeometry], maxsr, maxlr);
189 void init_neighbor_list(FILE *log, t_forcerec *fr, int homenr)
191 /* Make maxlr tunable! (does not seem to be a big difference though)
192 * This parameter determines the number of i particles in a long range
193 * neighbourlist. Too few means many function calls, too many means
196 int maxsr, maxsr_wat, maxlr, maxlr_wat;
197 int ielec, ielecf, ivdw, ielecmod, ielecmodf, ivdwmod, type;
199 int igeometry_def, igeometry_w, igeometry_ww;
203 /* maxsr = homenr-fr->nWatMol*3; */
208 gmx_fatal(FARGS, "%s, %d: Negative number of short range atoms.\n"
209 "Call your Gromacs dealer for assistance.", __FILE__, __LINE__);
211 /* This is just for initial allocation, so we do not reallocate
212 * all the nlist arrays many times in a row.
213 * The numbers seem very accurate, but they are uncritical.
215 maxsr_wat = min(fr->nWatMol, (homenr+2)/3);
219 maxlr_wat = min(maxsr_wat, maxlr);
223 maxlr = maxlr_wat = 0;
226 /* Determine the values for ielec/ivdw. */
227 ielec = fr->nbkernel_elec_interaction;
228 ivdw = fr->nbkernel_vdw_interaction;
229 ielecmod = fr->nbkernel_elec_modifier;
230 ivdwmod = fr->nbkernel_vdw_modifier;
231 type = GMX_NBLIST_INTERACTION_STANDARD;
233 fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
236 igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
240 igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
243 fprintf(log, "\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
247 if (fr->solvent_opt == esolTIP4P)
249 igeometry_w = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
250 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
254 igeometry_w = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
255 igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
258 for (i = 0; i < fr->nnblists; i++)
260 nbl = &(fr->nblists[i]);
262 if ((fr->adress_type != eAdressOff) && (i >= fr->nnblists/2))
264 type = GMX_NBLIST_INTERACTION_ADRESS;
266 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ], &nbl->nlist_lr[eNL_VDWQQ],
267 maxsr, maxlr, ivdw, ivdwmod, ielec, ielecmod, igeometry_def, type);
268 init_nblist(log, &nbl->nlist_sr[eNL_VDW], &nbl->nlist_lr[eNL_VDW],
269 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, igeometry_def, type);
270 init_nblist(log, &nbl->nlist_sr[eNL_QQ], &nbl->nlist_lr[eNL_QQ],
271 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_def, type);
272 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATER], &nbl->nlist_lr[eNL_VDWQQ_WATER],
273 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_w, type);
274 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATER], &nbl->nlist_lr[eNL_QQ_WATER],
275 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_w, type);
276 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_WATERWATER], &nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
277 maxsr_wat, maxlr_wat, ivdw, ivdwmod, ielec, ielecmod, igeometry_ww, type);
278 init_nblist(log, &nbl->nlist_sr[eNL_QQ_WATERWATER], &nbl->nlist_lr[eNL_QQ_WATERWATER],
279 maxsr_wat, maxlr_wat, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielec, ielecmod, igeometry_ww, type);
281 /* Did we get the solvent loops so we can use optimized water kernels? */
282 if (nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf == NULL
283 || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf == NULL
284 #ifndef DISABLE_WATERWATER_NLIST
285 || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf == NULL
286 || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf == NULL
290 fr->solvent_opt = esolNO;
291 fprintf(log, "Note: The available nonbonded kernels do not support water optimization - disabling.\n");
294 if (fr->efep != efepNO)
296 if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
298 ielecf = GMX_NBKERNEL_ELEC_EWALD;
299 ielecmodf = eintmodNONE;
304 ielecmodf = ielecmod;
307 init_nblist(log, &nbl->nlist_sr[eNL_VDWQQ_FREE], &nbl->nlist_lr[eNL_VDWQQ_FREE],
308 maxsr, maxlr, ivdw, ivdwmod, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
309 init_nblist(log, &nbl->nlist_sr[eNL_VDW_FREE], &nbl->nlist_lr[eNL_VDW_FREE],
310 maxsr, maxlr, ivdw, ivdwmod, GMX_NBKERNEL_ELEC_NONE, eintmodNONE, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
311 init_nblist(log, &nbl->nlist_sr[eNL_QQ_FREE], &nbl->nlist_lr[eNL_QQ_FREE],
312 maxsr, maxlr, GMX_NBKERNEL_VDW_NONE, eintmodNONE, ielecf, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
316 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
318 init_nblist(log, &fr->QMMMlist, NULL,
319 maxsr, maxlr, 0, 0, ielec, ielecmod, GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
327 fr->ns.nblist_initialized = TRUE;
330 static void reset_nblist(t_nblist *nl)
341 static void reset_neighbor_lists(t_forcerec *fr, gmx_bool bResetSR, gmx_bool bResetLR)
347 /* only reset the short-range nblist */
348 reset_nblist(&(fr->QMMMlist));
351 for (n = 0; n < fr->nnblists; n++)
353 for (i = 0; i < eNL_NR; i++)
357 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
361 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
370 static inline void new_i_nblist(t_nblist *nlist,
371 gmx_bool bLR, atom_id i_atom, int shift, int gid)
373 int i, k, nri, nshift;
377 /* Check whether we have to increase the i counter */
379 (nlist->iinr[nri] != i_atom) ||
380 (nlist->shift[nri] != shift) ||
381 (nlist->gid[nri] != gid))
383 /* This is something else. Now see if any entries have
384 * been added in the list of the previous atom.
387 ((nlist->jindex[nri+1] > nlist->jindex[nri]) &&
388 (nlist->gid[nri] != -1)))
390 /* If so increase the counter */
393 if (nlist->nri >= nlist->maxnri)
395 nlist->maxnri += over_alloc_large(nlist->nri);
396 reallocate_nblist(nlist);
399 /* Set the number of neighbours and the atom number */
400 nlist->jindex[nri+1] = nlist->jindex[nri];
401 nlist->iinr[nri] = i_atom;
402 nlist->gid[nri] = gid;
403 nlist->shift[nri] = shift;
407 static inline void close_i_nblist(t_nblist *nlist)
409 int nri = nlist->nri;
414 /* Add elements up to padding. Since we allocate memory in units
415 * of the simd_padding width, we do not have to check for possible
416 * list reallocation here.
418 while ((nlist->nrj % nlist->simd_padding_width) != 0)
420 /* Use -4 here, so we can write forces for 4 atoms before real data */
421 nlist->jjnr[nlist->nrj++] = -4;
423 nlist->jindex[nri+1] = nlist->nrj;
425 len = nlist->nrj - nlist->jindex[nri];
427 /* nlist length for water i molecules is treated statically
430 if (len > nlist->maxlen)
437 static inline void close_nblist(t_nblist *nlist)
439 /* Only close this nblist when it has been initialized.
440 * Avoid the creation of i-lists with no j-particles.
444 /* Some assembly kernels do not support empty lists,
445 * make sure here that we don't generate any empty lists.
446 * With the current ns code this branch is taken in two cases:
447 * No i-particles at all: nri=-1 here
448 * There are i-particles, but no j-particles; nri=0 here
454 /* Close list number nri by incrementing the count */
459 static inline void close_neighbor_lists(t_forcerec *fr, gmx_bool bMakeQMMMnblist)
465 close_nblist(&(fr->QMMMlist));
468 for (n = 0; n < fr->nnblists; n++)
470 for (i = 0; (i < eNL_NR); i++)
472 close_nblist(&(fr->nblists[n].nlist_sr[i]));
473 close_nblist(&(fr->nblists[n].nlist_lr[i]));
479 static inline void add_j_to_nblist(t_nblist *nlist, atom_id j_atom, gmx_bool bLR)
481 int nrj = nlist->nrj;
483 if (nlist->nrj >= nlist->maxnrj)
485 nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1), nlist->simd_padding_width);
489 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
490 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
493 srenew(nlist->jjnr, nlist->maxnrj);
496 nlist->jjnr[nrj] = j_atom;
500 static inline void add_j_to_nblist_cg(t_nblist *nlist,
501 atom_id j_start, int j_end,
502 t_excl *bexcl, gmx_bool i_is_j,
505 int nrj = nlist->nrj;
508 if (nlist->nrj >= nlist->maxnrj)
510 nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
513 fprintf(debug, "Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
514 bLR ? "LR" : "SR", nlist->ielec, nlist->ivdw, nlist->type, nlist->igeometry, nlist->maxnrj);
517 srenew(nlist->jjnr, nlist->maxnrj);
518 srenew(nlist->jjnr_end, nlist->maxnrj);
519 srenew(nlist->excl, nlist->maxnrj*MAX_CGCGSIZE);
522 nlist->jjnr[nrj] = j_start;
523 nlist->jjnr_end[nrj] = j_end;
525 if (j_end - j_start > MAX_CGCGSIZE)
527 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);
530 /* Set the exclusions */
531 for (j = j_start; j < j_end; j++)
533 nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
537 /* Avoid double counting of intra-cg interactions */
538 for (j = 1; j < j_end-j_start; j++)
540 nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
548 put_in_list_t (gmx_bool bHaveVdW[],
565 put_in_list_at(gmx_bool bHaveVdW[],
581 /* The a[] index has been removed,
582 * to put it back in i_atom should be a[i0] and jj should be a[jj].
587 t_nblist * vdwc_free = NULL;
588 t_nblist * vdw_free = NULL;
589 t_nblist * coul_free = NULL;
590 t_nblist * vdwc_ww = NULL;
591 t_nblist * coul_ww = NULL;
593 int i, j, jcg, igid, gid, nbl_ind, ind_ij;
594 atom_id jj, jj0, jj1, i_atom;
599 real *charge, *chargeB;
600 real qi, qiB, qq, rlj;
601 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
602 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
606 /* Copy some pointers */
608 charge = md->chargeA;
609 chargeB = md->chargeB;
612 bPert = md->bPerturbed;
616 nicg = index[icg+1]-i0;
618 /* Get the i charge group info */
619 igid = GET_CGINFO_GID(cginfo[icg]);
621 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
626 /* Check if any of the particles involved are perturbed.
627 * If not we can do the cheaper normal put_in_list
628 * and use more solvent optimization.
630 for (i = 0; i < nicg; i++)
632 bFreeEnergy |= bPert[i0+i];
634 /* Loop over the j charge groups */
635 for (j = 0; (j < nj && !bFreeEnergy); j++)
640 /* Finally loop over the atoms in the j-charge group */
641 for (jj = jj0; jj < jj1; jj++)
643 bFreeEnergy |= bPert[jj];
648 /* Unpack pointers to neighbourlist structs */
649 if (fr->nnblists == 1)
655 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
659 nlist = fr->nblists[nbl_ind].nlist_lr;
663 nlist = fr->nblists[nbl_ind].nlist_sr;
666 if (iwater != esolNO)
668 vdwc = &nlist[eNL_VDWQQ_WATER];
669 vdw = &nlist[eNL_VDW];
670 coul = &nlist[eNL_QQ_WATER];
671 #ifndef DISABLE_WATERWATER_NLIST
672 vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
673 coul_ww = &nlist[eNL_QQ_WATERWATER];
678 vdwc = &nlist[eNL_VDWQQ];
679 vdw = &nlist[eNL_VDW];
680 coul = &nlist[eNL_QQ];
685 if (iwater != esolNO)
687 /* Loop over the atoms in the i charge group */
689 gid = GID(igid, jgid, ngid);
690 /* Create new i_atom for each energy group */
691 if (bDoCoul && bDoVdW)
693 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
694 #ifndef DISABLE_WATERWATER_NLIST
695 new_i_nblist(vdwc_ww, bLR, i_atom, shift, gid);
700 new_i_nblist(vdw, bLR, i_atom, shift, gid);
704 new_i_nblist(coul, bLR, i_atom, shift, gid);
705 #ifndef DISABLE_WATERWATER_NLIST
706 new_i_nblist(coul_ww, bLR, i_atom, shift, gid);
709 /* Loop over the j charge groups */
710 for (j = 0; (j < nj); j++)
720 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
722 if (iwater == esolSPC && jwater == esolSPC)
724 /* Interaction between two SPC molecules */
727 /* VdW only - only first atoms in each water interact */
728 add_j_to_nblist(vdw, jj0, bLR);
732 #ifdef DISABLE_WATERWATER_NLIST
733 /* Add entries for the three atoms - only do VdW if we need to */
736 add_j_to_nblist(coul, jj0, bLR);
740 add_j_to_nblist(vdwc, jj0, bLR);
742 add_j_to_nblist(coul, jj0+1, bLR);
743 add_j_to_nblist(coul, jj0+2, bLR);
745 /* One entry for the entire water-water interaction */
748 add_j_to_nblist(coul_ww, jj0, bLR);
752 add_j_to_nblist(vdwc_ww, jj0, bLR);
757 else if (iwater == esolTIP4P && jwater == esolTIP4P)
759 /* Interaction between two TIP4p molecules */
762 /* VdW only - only first atoms in each water interact */
763 add_j_to_nblist(vdw, jj0, bLR);
767 #ifdef DISABLE_WATERWATER_NLIST
768 /* Add entries for the four atoms - only do VdW if we need to */
771 add_j_to_nblist(vdw, jj0, bLR);
773 add_j_to_nblist(coul, jj0+1, bLR);
774 add_j_to_nblist(coul, jj0+2, bLR);
775 add_j_to_nblist(coul, jj0+3, bLR);
777 /* One entry for the entire water-water interaction */
780 add_j_to_nblist(coul_ww, jj0, bLR);
784 add_j_to_nblist(vdwc_ww, jj0, bLR);
791 /* j charge group is not water, but i is.
792 * Add entries to the water-other_atom lists; the geometry of the water
793 * molecule doesn't matter - that is taken care of in the nonbonded kernel,
794 * so we don't care if it is SPC or TIP4P...
801 for (jj = jj0; (jj < jj1); jj++)
805 add_j_to_nblist(coul, jj, bLR);
811 for (jj = jj0; (jj < jj1); jj++)
813 if (bHaveVdW[type[jj]])
815 add_j_to_nblist(vdw, jj, bLR);
821 /* _charge_ _groups_ interact with both coulomb and LJ */
822 /* Check which atoms we should add to the lists! */
823 for (jj = jj0; (jj < jj1); jj++)
825 if (bHaveVdW[type[jj]])
829 add_j_to_nblist(vdwc, jj, bLR);
833 add_j_to_nblist(vdw, jj, bLR);
836 else if (charge[jj] != 0)
838 add_j_to_nblist(coul, jj, bLR);
845 close_i_nblist(coul);
846 close_i_nblist(vdwc);
847 #ifndef DISABLE_WATERWATER_NLIST
848 close_i_nblist(coul_ww);
849 close_i_nblist(vdwc_ww);
854 /* no solvent as i charge group */
855 /* Loop over the atoms in the i charge group */
856 for (i = 0; i < nicg; i++)
859 gid = GID(igid, jgid, ngid);
862 /* Create new i_atom for each energy group */
863 if (bDoVdW && bDoCoul)
865 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
869 new_i_nblist(vdw, bLR, i_atom, shift, gid);
873 new_i_nblist(coul, bLR, i_atom, shift, gid);
875 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
876 bDoCoul_i = (bDoCoul && qi != 0);
878 if (bDoVdW_i || bDoCoul_i)
880 /* Loop over the j charge groups */
881 for (j = 0; (j < nj); j++)
885 /* Check for large charge groups */
896 /* Finally loop over the atoms in the j-charge group */
897 for (jj = jj0; jj < jj1; jj++)
899 bNotEx = NOTEXCL(bExcl, i, jj);
907 add_j_to_nblist(coul, jj, bLR);
912 if (bHaveVdW[type[jj]])
914 add_j_to_nblist(vdw, jj, bLR);
919 if (bHaveVdW[type[jj]])
923 add_j_to_nblist(vdwc, jj, bLR);
927 add_j_to_nblist(vdw, jj, bLR);
930 else if (charge[jj] != 0)
932 add_j_to_nblist(coul, jj, bLR);
940 close_i_nblist(coul);
941 close_i_nblist(vdwc);
947 /* we are doing free energy */
948 vdwc_free = &nlist[eNL_VDWQQ_FREE];
949 vdw_free = &nlist[eNL_VDW_FREE];
950 coul_free = &nlist[eNL_QQ_FREE];
951 /* Loop over the atoms in the i charge group */
952 for (i = 0; i < nicg; i++)
955 gid = GID(igid, jgid, ngid);
957 qiB = chargeB[i_atom];
959 /* Create new i_atom for each energy group */
960 if (bDoVdW && bDoCoul)
962 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
966 new_i_nblist(vdw, bLR, i_atom, shift, gid);
970 new_i_nblist(coul, bLR, i_atom, shift, gid);
973 new_i_nblist(vdw_free, bLR, i_atom, shift, gid);
974 new_i_nblist(coul_free, bLR, i_atom, shift, gid);
975 new_i_nblist(vdwc_free, bLR, i_atom, shift, gid);
977 bDoVdW_i = (bDoVdW &&
978 (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
979 bDoCoul_i = (bDoCoul && (qi != 0 || qiB != 0));
980 /* For TIP4P the first atom does not have a charge,
981 * but the last three do. So we should still put an atom
982 * without LJ but with charge in the water-atom neighborlist
983 * for a TIP4p i charge group.
984 * For SPC type water the first atom has LJ and charge,
985 * so there is no such problem.
987 if (iwater == esolNO)
989 bDoCoul_i_sol = bDoCoul_i;
993 bDoCoul_i_sol = bDoCoul;
996 if (bDoVdW_i || bDoCoul_i_sol)
998 /* Loop over the j charge groups */
999 for (j = 0; (j < nj); j++)
1003 /* Check for large charge groups */
1014 /* Finally loop over the atoms in the j-charge group */
1015 bFree = bPert[i_atom];
1016 for (jj = jj0; (jj < jj1); jj++)
1018 bFreeJ = bFree || bPert[jj];
1019 /* Complicated if, because the water H's should also
1020 * see perturbed j-particles
1022 if (iwater == esolNO || i == 0 || bFreeJ)
1024 bNotEx = NOTEXCL(bExcl, i, jj);
1032 if (charge[jj] != 0 || chargeB[jj] != 0)
1034 add_j_to_nblist(coul_free, jj, bLR);
1037 else if (!bDoCoul_i)
1039 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1041 add_j_to_nblist(vdw_free, jj, bLR);
1046 if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1048 if (charge[jj] != 0 || chargeB[jj] != 0)
1050 add_j_to_nblist(vdwc_free, jj, bLR);
1054 add_j_to_nblist(vdw_free, jj, bLR);
1057 else if (charge[jj] != 0 || chargeB[jj] != 0)
1059 add_j_to_nblist(coul_free, jj, bLR);
1065 /* This is done whether or not bWater is set */
1066 if (charge[jj] != 0)
1068 add_j_to_nblist(coul, jj, bLR);
1071 else if (!bDoCoul_i_sol)
1073 if (bHaveVdW[type[jj]])
1075 add_j_to_nblist(vdw, jj, bLR);
1080 if (bHaveVdW[type[jj]])
1082 if (charge[jj] != 0)
1084 add_j_to_nblist(vdwc, jj, bLR);
1088 add_j_to_nblist(vdw, jj, bLR);
1091 else if (charge[jj] != 0)
1093 add_j_to_nblist(coul, jj, bLR);
1101 close_i_nblist(vdw);
1102 close_i_nblist(coul);
1103 close_i_nblist(vdwc);
1104 close_i_nblist(vdw_free);
1105 close_i_nblist(coul_free);
1106 close_i_nblist(vdwc_free);
1112 put_in_list_adress(gmx_bool bHaveVdW[],
1128 /* The a[] index has been removed,
1129 * to put it back in i_atom should be a[i0] and jj should be a[jj].
1134 t_nblist * vdwc_adress = NULL;
1135 t_nblist * vdw_adress = NULL;
1136 t_nblist * coul_adress = NULL;
1137 t_nblist * vdwc_ww = NULL;
1138 t_nblist * coul_ww = NULL;
1140 int i, j, jcg, igid, gid, nbl_ind, nbl_ind_adress;
1141 atom_id jj, jj0, jj1, i_atom;
1146 real *charge, *chargeB;
1148 real qi, qiB, qq, rlj;
1149 gmx_bool bFreeEnergy, bFree, bFreeJ, bNotEx, *bPert;
1150 gmx_bool bDoVdW_i, bDoCoul_i, bDoCoul_i_sol;
1152 gmx_bool j_all_atom;
1154 t_nblist *nlist, *nlist_adress;
1156 /* Copy some pointers */
1157 cginfo = fr->cginfo;
1158 charge = md->chargeA;
1159 chargeB = md->chargeB;
1162 bPert = md->bPerturbed;
1165 /* Get atom range */
1167 nicg = index[icg+1]-i0;
1169 /* Get the i charge group info */
1170 igid = GET_CGINFO_GID(cginfo[icg]);
1172 iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1176 gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1179 /* Unpack pointers to neighbourlist structs */
1180 if (fr->nnblists == 2)
1187 nbl_ind = fr->gid2nblists[GID(igid, jgid, ngid)];
1188 nbl_ind_adress = nbl_ind+fr->nnblists/2;
1192 nlist = fr->nblists[nbl_ind].nlist_lr;
1193 nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1197 nlist = fr->nblists[nbl_ind].nlist_sr;
1198 nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1202 vdwc = &nlist[eNL_VDWQQ];
1203 vdw = &nlist[eNL_VDW];
1204 coul = &nlist[eNL_QQ];
1206 vdwc_adress = &nlist_adress[eNL_VDWQQ];
1207 vdw_adress = &nlist_adress[eNL_VDW];
1208 coul_adress = &nlist_adress[eNL_QQ];
1210 /* We do not support solvent optimization with AdResS for now.
1211 For this we would need hybrid solvent-other kernels */
1213 /* no solvent as i charge group */
1214 /* Loop over the atoms in the i charge group */
1215 for (i = 0; i < nicg; i++)
1218 gid = GID(igid, jgid, ngid);
1219 qi = charge[i_atom];
1221 /* Create new i_atom for each energy group */
1222 if (bDoVdW && bDoCoul)
1224 new_i_nblist(vdwc, bLR, i_atom, shift, gid);
1225 new_i_nblist(vdwc_adress, bLR, i_atom, shift, gid);
1230 new_i_nblist(vdw, bLR, i_atom, shift, gid);
1231 new_i_nblist(vdw_adress, bLR, i_atom, shift, gid);
1236 new_i_nblist(coul, bLR, i_atom, shift, gid);
1237 new_i_nblist(coul_adress, bLR, i_atom, shift, gid);
1239 bDoVdW_i = (bDoVdW && bHaveVdW[type[i_atom]]);
1240 bDoCoul_i = (bDoCoul && qi != 0);
1242 if (bDoVdW_i || bDoCoul_i)
1244 /* Loop over the j charge groups */
1245 for (j = 0; (j < nj); j++)
1249 /* Check for large charge groups */
1260 /* Finally loop over the atoms in the j-charge group */
1261 for (jj = jj0; jj < jj1; jj++)
1263 bNotEx = NOTEXCL(bExcl, i, jj);
1265 b_hybrid = !((wf[i_atom] == 1 && wf[jj] == 1) || (wf[i_atom] == 0 && wf[jj] == 0));
1271 if (charge[jj] != 0)
1275 add_j_to_nblist(coul, jj, bLR);
1279 add_j_to_nblist(coul_adress, jj, bLR);
1283 else if (!bDoCoul_i)
1285 if (bHaveVdW[type[jj]])
1289 add_j_to_nblist(vdw, jj, bLR);
1293 add_j_to_nblist(vdw_adress, jj, bLR);
1299 if (bHaveVdW[type[jj]])
1301 if (charge[jj] != 0)
1305 add_j_to_nblist(vdwc, jj, bLR);
1309 add_j_to_nblist(vdwc_adress, jj, bLR);
1316 add_j_to_nblist(vdw, jj, bLR);
1320 add_j_to_nblist(vdw_adress, jj, bLR);
1325 else if (charge[jj] != 0)
1329 add_j_to_nblist(coul, jj, bLR);
1333 add_j_to_nblist(coul_adress, jj, bLR);
1342 close_i_nblist(vdw);
1343 close_i_nblist(coul);
1344 close_i_nblist(vdwc);
1345 close_i_nblist(vdw_adress);
1346 close_i_nblist(coul_adress);
1347 close_i_nblist(vdwc_adress);
1353 put_in_list_qmmm(gmx_bool bHaveVdW[],
1370 int i, j, jcg, igid, gid;
1371 atom_id jj, jj0, jj1, i_atom;
1375 /* Get atom range */
1377 nicg = index[icg+1]-i0;
1379 /* Get the i charge group info */
1380 igid = GET_CGINFO_GID(fr->cginfo[icg]);
1382 coul = &fr->QMMMlist;
1384 /* Loop over atoms in the ith charge group */
1385 for (i = 0; i < nicg; i++)
1388 gid = GID(igid, jgid, ngid);
1389 /* Create new i_atom for each energy group */
1390 new_i_nblist(coul, bLR, i_atom, shift, gid);
1392 /* Loop over the j charge groups */
1393 for (j = 0; j < nj; j++)
1397 /* Charge groups cannot have QM and MM atoms simultaneously */
1402 /* Finally loop over the atoms in the j-charge group */
1403 for (jj = jj0; jj < jj1; jj++)
1405 bNotEx = NOTEXCL(bExcl, i, jj);
1408 add_j_to_nblist(coul, jj, bLR);
1413 close_i_nblist(coul);
1418 put_in_list_cg(gmx_bool bHaveVdW[],
1435 int igid, gid, nbl_ind;
1439 cginfo = fr->cginfo[icg];
1441 igid = GET_CGINFO_GID(cginfo);
1442 gid = GID(igid, jgid, ngid);
1444 /* Unpack pointers to neighbourlist structs */
1445 if (fr->nnblists == 1)
1451 nbl_ind = fr->gid2nblists[gid];
1455 vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1459 vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1462 /* Make a new neighbor list for charge group icg.
1463 * Currently simply one neighbor list is made with LJ and Coulomb.
1464 * If required, zero interactions could be removed here
1465 * or in the force loop.
1467 new_i_nblist(vdwc, bLR, index[icg], shift, gid);
1468 vdwc->iinr_end[vdwc->nri] = index[icg+1];
1470 for (j = 0; (j < nj); j++)
1473 /* Skip the icg-icg pairs if all self interactions are excluded */
1474 if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1476 /* Here we add the j charge group jcg to the list,
1477 * exclusions are also added to the list.
1479 add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1483 close_i_nblist(vdwc);
1486 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1493 for (i = start; i < end; i++)
1495 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1497 SETEXCL(bexcl, i-start, excl->a[k]);
1503 for (i = start; i < end; i++)
1505 for (k = excl->index[i]; k < excl->index[i+1]; k++)
1507 RMEXCL(bexcl, i-start, excl->a[k]);
1513 int calc_naaj(int icg, int cgtot)
1517 if ((cgtot % 2) == 1)
1519 /* Odd number of charge groups, easy */
1520 naaj = 1 + (cgtot/2);
1522 else if ((cgtot % 4) == 0)
1524 /* Multiple of four is hard */
1561 fprintf(log, "naaj=%d\n", naaj);
1567 /************************************************
1569 * S I M P L E C O R E S T U F F
1571 ************************************************/
1573 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1574 rvec b_inv, int *shift)
1576 /* This code assumes that the cut-off is smaller than
1577 * a half times the smallest diagonal element of the box.
1579 const real h25 = 2.5;
1584 /* Compute diff vector */
1585 dz = xj[ZZ] - xi[ZZ];
1586 dy = xj[YY] - xi[YY];
1587 dx = xj[XX] - xi[XX];
1589 /* Perform NINT operation, using trunc operation, therefore
1590 * we first add 2.5 then subtract 2 again
1592 tz = dz*b_inv[ZZ] + h25;
1594 dz -= tz*box[ZZ][ZZ];
1595 dy -= tz*box[ZZ][YY];
1596 dx -= tz*box[ZZ][XX];
1598 ty = dy*b_inv[YY] + h25;
1600 dy -= ty*box[YY][YY];
1601 dx -= ty*box[YY][XX];
1603 tx = dx*b_inv[XX]+h25;
1605 dx -= tx*box[XX][XX];
1607 /* Distance squared */
1608 r2 = (dx*dx) + (dy*dy) + (dz*dz);
1610 *shift = XYZ2IS(tx, ty, tz);
1615 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1616 rvec b_inv, int *shift)
1618 const real h15 = 1.5;
1624 /* Compute diff vector */
1625 dx = xj[XX] - xi[XX];
1626 dy = xj[YY] - xi[YY];
1627 dz = xj[ZZ] - xi[ZZ];
1629 /* Perform NINT operation, using trunc operation, therefore
1630 * we first add 1.5 then subtract 1 again
1632 tx = dx*b_inv[XX] + h15;
1633 ty = dy*b_inv[YY] + h15;
1634 tz = dz*b_inv[ZZ] + h15;
1639 /* Correct diff vector for translation */
1640 ddx = tx*box_size[XX] - dx;
1641 ddy = ty*box_size[YY] - dy;
1642 ddz = tz*box_size[ZZ] - dz;
1644 /* Distance squared */
1645 r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1647 *shift = XYZ2IS(tx, ty, tz);
1652 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1653 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1654 int icg, int jgid, t_block *cgs, t_excl bexcl[],
1655 int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1657 if (nsbuf->nj + nrj > MAX_CG)
1659 put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1660 cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1661 /* Reset buffer contents */
1662 nsbuf->ncg = nsbuf->nj = 0;
1664 nsbuf->jcg[nsbuf->ncg++] = cg_j;
1668 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1669 int njcg, atom_id jcg[],
1670 matrix box, rvec b_inv, real rcut2,
1671 t_block *cgs, t_ns_buf **ns_buf,
1672 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1673 t_excl bexcl[], t_forcerec *fr,
1674 put_in_list_t *put_in_list)
1678 int *cginfo = fr->cginfo;
1679 atom_id cg_j, *cgindex;
1682 cgindex = cgs->index;
1684 for (j = 0; (j < njcg); j++)
1687 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1688 if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1690 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1691 if (!(i_egp_flags[jgid] & EGP_EXCL))
1693 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1694 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1701 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1702 int njcg, atom_id jcg[],
1703 gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1704 t_block *cgs, t_ns_buf **ns_buf,
1705 gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1706 t_excl bexcl[], t_forcerec *fr,
1707 put_in_list_t *put_in_list)
1711 int *cginfo = fr->cginfo;
1712 atom_id cg_j, *cgindex;
1715 cgindex = cgs->index;
1719 for (j = 0; (j < njcg); j++)
1722 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1723 if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1725 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1726 if (!(i_egp_flags[jgid] & EGP_EXCL))
1728 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1729 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1737 for (j = 0; (j < njcg); j++)
1740 nrj = cgindex[cg_j+1]-cgindex[cg_j];
1741 if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1743 jgid = GET_CGINFO_GID(cginfo[cg_j]);
1744 if (!(i_egp_flags[jgid] & EGP_EXCL))
1746 add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1747 bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1755 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1757 static int ns_simple_core(t_forcerec *fr,
1758 gmx_localtop_t *top,
1760 matrix box, rvec box_size,
1761 t_excl bexcl[], atom_id *aaj,
1762 int ngid, t_ns_buf **ns_buf,
1763 put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1767 int nsearch, icg, jcg, igid, i0, nri, nn;
1770 /* atom_id *i_atoms; */
1771 t_block *cgs = &(top->cgs);
1772 t_blocka *excl = &(top->excls);
1775 gmx_bool bBox, bTriclinic;
1778 rlist2 = sqr(fr->rlist);
1780 bBox = (fr->ePBC != epbcNONE);
1783 for (m = 0; (m < DIM); m++)
1785 b_inv[m] = divide_err(1.0, box_size[m]);
1787 bTriclinic = TRICLINIC(box);
1794 cginfo = fr->cginfo;
1797 for (icg = fr->cg0; (icg < fr->hcg); icg++)
1800 i0 = cgs->index[icg];
1801 nri = cgs->index[icg+1]-i0;
1802 i_atoms = &(cgs->a[i0]);
1803 i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1804 setexcl(nri,i_atoms,excl,TRUE,bexcl);
1806 igid = GET_CGINFO_GID(cginfo[icg]);
1807 i_egp_flags = fr->egp_flags + ngid*igid;
1808 setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1810 naaj = calc_naaj(icg, cgs->nr);
1813 ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1814 box, b_inv, rlist2, cgs, ns_buf,
1815 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1819 ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1820 bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1821 bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1825 for (nn = 0; (nn < ngid); nn++)
1827 for (k = 0; (k < SHIFTS); k++)
1829 nsbuf = &(ns_buf[nn][k]);
1832 put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1833 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1834 nsbuf->ncg = nsbuf->nj = 0;
1838 /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1839 setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1841 close_neighbor_lists(fr, FALSE);
1846 /************************************************
1848 * N S 5 G R I D S T U F F
1850 ************************************************/
1852 static inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1853 int *dx0, int *dx1, real *dcx2)
1881 for (i = xgi0; i >= 0; i--)
1883 dcx = (i+1)*gridx-x;
1892 for (i = xgi1; i < Nx; i++)
1905 static inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1906 int ncpddc, int shift_min, int shift_max,
1907 int *g0, int *g1, real *dcx2)
1910 int g_min, g_max, shift_home;
1943 g_min = (shift_min == shift_home ? 0 : ncpddc);
1944 g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1951 else if (shift_max < 0)
1966 /* Check one grid cell down */
1967 dcx = ((*g0 - 1) + 1)*gridx - x;
1979 /* Check one grid cell up */
1980 dcx = (*g1 + 1)*gridx - x;
1992 #define sqr(x) ((x)*(x))
1993 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1994 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1995 /****************************************************
1997 * F A S T N E I G H B O R S E A R C H I N G
1999 * Optimized neighboursearching routine using grid
2000 * at least 1x1x1, see GROMACS manual
2002 ****************************************************/
2005 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2006 real *rvdw2, real *rcoul2,
2007 real *rs2, real *rm2, real *rl2)
2009 *rs2 = sqr(fr->rlist);
2011 if (bDoLongRange && fr->bTwinRange)
2013 /* The VdW and elec. LR cut-off's could be different,
2014 * so we can not simply set them to rlistlong.
2016 if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
2017 fr->rvdw > fr->rlist)
2019 *rvdw2 = sqr(fr->rlistlong);
2023 *rvdw2 = sqr(fr->rvdw);
2025 if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
2026 fr->rcoulomb > fr->rlist)
2028 *rcoul2 = sqr(fr->rlistlong);
2032 *rcoul2 = sqr(fr->rcoulomb);
2037 /* Workaround for a gcc -O3 or -ffast-math problem */
2041 *rm2 = min(*rvdw2, *rcoul2);
2042 *rl2 = max(*rvdw2, *rcoul2);
2045 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2047 real rvdw2, rcoul2, rs2, rm2, rl2;
2050 get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2052 /* Short range buffers */
2053 snew(ns->nl_sr, ngid);
2055 snew(ns->nsr, ngid);
2056 snew(ns->nlr_ljc, ngid);
2057 snew(ns->nlr_one, ngid);
2059 /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2060 /* Long range VdW and Coul buffers */
2061 snew(ns->nl_lr_ljc, ngid);
2062 /* Long range VdW or Coul only buffers */
2063 snew(ns->nl_lr_one, ngid);
2065 for (j = 0; (j < ngid); j++)
2067 snew(ns->nl_sr[j], MAX_CG);
2068 snew(ns->nl_lr_ljc[j], MAX_CG);
2069 snew(ns->nl_lr_one[j], MAX_CG);
2074 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2079 static int nsgrid_core(FILE *log, t_commrec *cr, t_forcerec *fr,
2080 matrix box, rvec box_size, int ngid,
2081 gmx_localtop_t *top,
2082 t_grid *grid, rvec x[],
2083 t_excl bexcl[], gmx_bool *bExcludeAlleg,
2084 t_nrnb *nrnb, t_mdatoms *md,
2085 real *lambda, real *dvdlambda,
2086 gmx_grppairener_t *grppener,
2087 put_in_list_t *put_in_list,
2088 gmx_bool bHaveVdW[],
2089 gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2092 atom_id **nl_lr_ljc, **nl_lr_one, **nl_sr;
2093 int *nlr_ljc, *nlr_one, *nsr;
2094 gmx_domdec_t *dd = NULL;
2095 t_block *cgs = &(top->cgs);
2096 int *cginfo = fr->cginfo;
2097 /* atom_id *i_atoms,*cgsindex=cgs->index; */
2099 int cell_x, cell_y, cell_z;
2100 int d, tx, ty, tz, dx, dy, dz, cj;
2101 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2102 int zsh_ty, zsh_tx, ysh_tx;
2104 int dx0, dx1, dy0, dy1, dz0, dz1;
2105 int Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2106 real gridx, gridy, gridz, grid_x, grid_y, grid_z;
2107 real *dcx2, *dcy2, *dcz2;
2109 int cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2110 int jcg0, jcg1, jjcg, cgj0, jgid;
2111 int *grida, *gridnra, *gridind;
2112 gmx_bool rvdw_lt_rcoul, rcoul_lt_rvdw;
2113 rvec xi, *cgcm, grid_offset;
2114 real r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2116 gmx_bool bDomDec, bTriclinicX, bTriclinicY;
2121 bDomDec = DOMAINDECOMP(cr);
2127 bTriclinicX = ((YY < grid->npbcdim &&
2128 (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2129 (ZZ < grid->npbcdim &&
2130 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2131 bTriclinicY = (ZZ < grid->npbcdim &&
2132 (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2136 get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2138 rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2139 rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2141 if (bMakeQMMMnblist)
2149 nl_lr_ljc = ns->nl_lr_ljc;
2150 nl_lr_one = ns->nl_lr_one;
2151 nlr_ljc = ns->nlr_ljc;
2152 nlr_one = ns->nlr_one;
2160 gridind = grid->index;
2161 gridnra = grid->nra;
2164 gridx = grid->cell_size[XX];
2165 gridy = grid->cell_size[YY];
2166 gridz = grid->cell_size[ZZ];
2170 copy_rvec(grid->cell_offset, grid_offset);
2171 copy_ivec(grid->ncpddc, ncpddc);
2176 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2177 zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2178 zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2179 ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2180 if (zsh_tx != 0 && ysh_tx != 0)
2182 /* This could happen due to rounding, when both ratios are 0.5 */
2191 /* We only want a list for the test particle */
2200 /* Set the shift range */
2201 for (d = 0; d < DIM; d++)
2205 /* Check if we need periodicity shifts.
2206 * Without PBC or with domain decomposition we don't need them.
2208 if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2215 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2226 /* Loop over charge groups */
2227 for (icg = cg0; (icg < cg1); icg++)
2229 igid = GET_CGINFO_GID(cginfo[icg]);
2230 /* Skip this charge group if all energy groups are excluded! */
2231 if (bExcludeAlleg[igid])
2236 i0 = cgs->index[icg];
2238 if (bMakeQMMMnblist)
2240 /* Skip this charge group if it is not a QM atom while making a
2241 * QM/MM neighbourlist
2243 if (md->bQM[i0] == FALSE)
2245 continue; /* MM particle, go to next particle */
2248 /* Compute the number of charge groups that fall within the control
2251 naaj = calc_naaj(icg, cgsnr);
2258 /* make a normal neighbourlist */
2262 /* Get the j charge-group and dd cell shift ranges */
2263 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2268 /* Compute the number of charge groups that fall within the control
2271 naaj = calc_naaj(icg, cgsnr);
2277 /* The i-particle is awlways the test particle,
2278 * so we want all j-particles
2280 max_jcg = cgsnr - 1;
2284 max_jcg = jcg1 - cgsnr;
2289 i_egp_flags = fr->egp_flags + igid*ngid;
2291 /* Set the exclusions for the atoms in charge group icg using a bitmask */
2292 setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2294 ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2296 /* Changed iicg to icg, DvdS 990115
2297 * (but see consistency check above, DvdS 990330)
2300 fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2301 icg, naaj, cell_x, cell_y, cell_z);
2303 /* Loop over shift vectors in three dimensions */
2304 for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2306 ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2307 /* Calculate range of cells in Z direction that have the shift tz */
2308 zgi = cell_z + tz*Nz;
2311 get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2313 get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2314 ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2320 for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2322 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2323 /* Calculate range of cells in Y direction that have the shift ty */
2326 ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2330 ygi = cell_y + ty*Ny;
2333 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2335 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2336 ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2342 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2344 XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2345 /* Calculate range of cells in X direction that have the shift tx */
2348 xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2352 xgi = cell_x + tx*Nx;
2355 get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2357 get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2358 ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2364 /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2365 * from the neigbour list as it will not interact */
2366 if (fr->adress_type != eAdressOff)
2368 if (md->wf[cgs->index[icg]] == 0 && egp_explicit(fr, igid))
2373 /* Get shift vector */
2374 shift = XYZ2IS(tx, ty, tz);
2376 range_check(shift, 0, SHIFTS);
2378 for (nn = 0; (nn < ngid); nn++)
2385 fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2386 shift, dx0, dx1, dy0, dy1, dz0, dz1);
2387 fprintf(log, "cgcm: %8.3f %8.3f %8.3f\n", cgcm[icg][XX],
2388 cgcm[icg][YY], cgcm[icg][ZZ]);
2389 fprintf(log, "xi: %8.3f %8.3f %8.3f\n", XI, YI, ZI);
2391 for (dx = dx0; (dx <= dx1); dx++)
2393 tmp1 = rl2 - dcx2[dx];
2394 for (dy = dy0; (dy <= dy1); dy++)
2396 tmp2 = tmp1 - dcy2[dy];
2399 for (dz = dz0; (dz <= dz1); dz++)
2401 if (tmp2 > dcz2[dz])
2403 /* Find grid-cell cj in which possible neighbours are */
2404 cj = xyz2ci(Ny, Nz, dx, dy, dz);
2406 /* Check out how many cgs (nrj) there in this cell */
2409 /* Find the offset in the cg list */
2412 /* Check if all j's are out of range so we
2413 * can skip the whole cell.
2414 * Should save some time, especially with DD.
2417 (grida[cgj0] >= max_jcg &&
2418 (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2424 for (j = 0; (j < nrj); j++)
2426 jjcg = grida[cgj0+j];
2428 /* check whether this guy is in range! */
2429 if ((jjcg >= jcg0 && jjcg < jcg1) ||
2432 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2435 /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2436 jgid = GET_CGINFO_GID(cginfo[jjcg]);
2437 /* check energy group exclusions */
2438 if (!(i_egp_flags[jgid] & EGP_EXCL))
2442 if (nsr[jgid] >= MAX_CG)
2444 /* Add to short-range list */
2445 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2446 nsr[jgid], nl_sr[jgid],
2447 cgs->index, /* cgsatoms, */ bexcl,
2448 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2451 nl_sr[jgid][nsr[jgid]++] = jjcg;
2455 if (nlr_ljc[jgid] >= MAX_CG)
2457 /* Add to LJ+coulomb long-range list */
2458 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2459 nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2460 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2463 nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2467 if (nlr_one[jgid] >= MAX_CG)
2469 /* Add to long-range list with only coul, or only LJ */
2470 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2471 nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2472 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2475 nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2487 /* CHECK whether there is anything left in the buffers */
2488 for (nn = 0; (nn < ngid); nn++)
2492 put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2493 cgs->index, /* cgsatoms, */ bexcl,
2494 shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2497 if (nlr_ljc[nn] > 0)
2499 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2500 nl_lr_ljc[nn], top->cgs.index,
2501 bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2504 if (nlr_one[nn] > 0)
2506 put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2507 nl_lr_one[nn], top->cgs.index,
2508 bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2514 /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2515 setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2517 /* No need to perform any left-over force calculations anymore (as we used to do here)
2518 * since we now save the proper long-range lists for later evaluation.
2523 /* Close neighbourlists */
2524 close_neighbor_lists(fr, bMakeQMMMnblist);
2529 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2533 if (natoms > ns->nra_alloc)
2535 ns->nra_alloc = over_alloc_dd(natoms);
2536 srenew(ns->bexcl, ns->nra_alloc);
2537 for (i = 0; i < ns->nra_alloc; i++)
2544 void init_ns(FILE *fplog, const t_commrec *cr,
2545 gmx_ns_t *ns, t_forcerec *fr,
2546 const gmx_mtop_t *mtop,
2549 int mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2553 /* Compute largest charge groups size (# atoms) */
2555 for (mt = 0; mt < mtop->nmoltype; mt++)
2557 cgs = &mtop->moltype[mt].cgs;
2558 for (icg = 0; (icg < cgs->nr); icg++)
2560 nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2564 /* Verify whether largest charge group is <= max cg.
2565 * This is determined by the type of the local exclusion type
2566 * Exclusions are stored in bits. (If the type is not large
2567 * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2569 maxcg = sizeof(t_excl)*8;
2570 if (nr_in_cg > maxcg)
2572 gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2576 ngid = mtop->groups.grps[egcENER].nr;
2577 snew(ns->bExcludeAlleg, ngid);
2578 for (i = 0; i < ngid; i++)
2580 ns->bExcludeAlleg[i] = TRUE;
2581 for (j = 0; j < ngid; j++)
2583 if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2585 ns->bExcludeAlleg[i] = FALSE;
2593 ns->grid = init_grid(fplog, fr);
2594 init_nsgrid_lists(fr, ngid, ns);
2599 snew(ns->ns_buf, ngid);
2600 for (i = 0; (i < ngid); i++)
2602 snew(ns->ns_buf[i], SHIFTS);
2604 ncg = ncg_mtop(mtop);
2605 snew(ns->simple_aaj, 2*ncg);
2606 for (jcg = 0; (jcg < ncg); jcg++)
2608 ns->simple_aaj[jcg] = jcg;
2609 ns->simple_aaj[jcg+ncg] = jcg;
2613 /* Create array that determines whether or not atoms have VdW */
2614 snew(ns->bHaveVdW, fr->ntype);
2615 for (i = 0; (i < fr->ntype); i++)
2617 for (j = 0; (j < fr->ntype); j++)
2619 ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2621 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2622 (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2623 (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2624 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2625 (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2630 pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2635 if (!DOMAINDECOMP(cr))
2637 /* This could be reduced with particle decomposition */
2638 ns_realloc_natoms(ns, mtop->natoms);
2641 ns->nblist_initialized = FALSE;
2643 /* nbr list debug dump */
2645 char *ptr = getenv("GMX_DUMP_NL");
2648 ns->dump_nl = strtol(ptr, NULL, 10);
2651 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2662 int search_neighbours(FILE *log, t_forcerec *fr,
2663 rvec x[], matrix box,
2664 gmx_localtop_t *top,
2665 gmx_groups_t *groups,
2667 t_nrnb *nrnb, t_mdatoms *md,
2668 real *lambda, real *dvdlambda,
2669 gmx_grppairener_t *grppener,
2671 gmx_bool bDoLongRangeNS,
2672 gmx_bool bPadListsForKernels)
2674 t_block *cgs = &(top->cgs);
2675 rvec box_size, grid_x0, grid_x1;
2677 real min_size, grid_dens;
2681 gmx_bool *i_egp_flags;
2682 int cg_start, cg_end, start, end;
2685 gmx_domdec_zones_t *dd_zones;
2686 put_in_list_t *put_in_list;
2690 /* Set some local variables */
2692 ngid = groups->grps[egcENER].nr;
2694 for (m = 0; (m < DIM); m++)
2696 box_size[m] = box[m][m];
2699 if (fr->ePBC != epbcNONE)
2701 if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2703 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.");
2707 min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2708 if (2*fr->rlistlong >= min_size)
2710 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2715 if (DOMAINDECOMP(cr))
2717 ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2721 /* Reset the neighbourlists */
2722 reset_neighbor_lists(fr, TRUE, TRUE);
2724 if (bGrid && bFillGrid)
2728 if (DOMAINDECOMP(cr))
2730 dd_zones = domdec_zones(cr->dd);
2736 get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2737 cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2739 grid_first(log, grid, NULL, NULL, fr->ePBC, box, grid_x0, grid_x1,
2740 fr->rlistlong, grid_dens);
2744 /* Don't know why this all is... (DvdS 3/99) */
2750 end = (cgs->nr+1)/2;
2753 if (DOMAINDECOMP(cr))
2756 fill_grid(log, dd_zones, grid, end, -1, end, fr->cg_cm);
2758 grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2762 fill_grid(log, NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2763 grid->icg0 = fr->cg0;
2764 grid->icg1 = fr->hcg;
2774 calc_elemnr(log, grid, start, end, cgs->nr);
2776 grid_last(log, grid, start, end, cgs->nr);
2780 check_grid(debug, grid);
2781 print_grid(debug, grid);
2786 /* Set the grid cell index for the test particle only.
2787 * The cell to cg index is not corrected, but that does not matter.
2789 fill_grid(log, NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2793 if (fr->adress_type == eAdressOff)
2795 if (!fr->ns.bCGlist)
2797 put_in_list = put_in_list_at;
2801 put_in_list = put_in_list_cg;
2806 put_in_list = put_in_list_adress;
2813 nsearch = nsgrid_core(log, cr, fr, box, box_size, ngid, top,
2814 grid, x, ns->bexcl, ns->bExcludeAlleg,
2815 nrnb, md, lambda, dvdlambda, grppener,
2816 put_in_list, ns->bHaveVdW,
2817 bDoLongRangeNS, FALSE);
2819 /* neighbour searching withouth QMMM! QM atoms have zero charge in
2820 * the classical calculation. The charge-charge interaction
2821 * between QM and MM atoms is handled in the QMMM core calculation
2822 * (see QMMM.c). The VDW however, we'd like to compute classically
2823 * and the QM MM atom pairs have just been put in the
2824 * corresponding neighbourlists. in case of QMMM we still need to
2825 * fill a special QMMM neighbourlist that contains all neighbours
2826 * of the QM atoms. If bQMMM is true, this list will now be made:
2828 if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2830 nsearch += nsgrid_core(log, cr, fr, box, box_size, ngid, top,
2831 grid, x, ns->bexcl, ns->bExcludeAlleg,
2832 nrnb, md, lambda, dvdlambda, grppener,
2833 put_in_list_qmmm, ns->bHaveVdW,
2834 bDoLongRangeNS, TRUE);
2839 nsearch = ns_simple_core(fr, top, md, box, box_size,
2840 ns->bexcl, ns->simple_aaj,
2841 ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2849 inc_nrnb(nrnb, eNR_NS, nsearch);
2850 /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2855 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2856 matrix scale_tot, rvec *x)
2858 int cg0, cg1, cg, a0, a1, a, i, j;
2859 real rint, hbuf2, scale;
2861 gmx_bool bIsotropic;
2866 rint = max(ir->rcoulomb, ir->rvdw);
2867 if (ir->rlist < rint)
2869 gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2877 if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2879 hbuf2 = sqr(0.5*(ir->rlist - rint));
2880 for (cg = cg0; cg < cg1; cg++)
2882 a0 = cgs->index[cg];
2883 a1 = cgs->index[cg+1];
2884 for (a = a0; a < a1; a++)
2886 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2896 scale = scale_tot[0][0];
2897 for (i = 1; i < DIM; i++)
2899 /* With anisotropic scaling, the original spherical ns volumes become
2900 * ellipsoids. To avoid costly transformations we use the minimum
2901 * eigenvalue of the scaling matrix for determining the buffer size.
2902 * Since the lower half is 0, the eigenvalues are the diagonal elements.
2904 scale = min(scale, scale_tot[i][i]);
2905 if (scale_tot[i][i] != scale_tot[i-1][i-1])
2909 for (j = 0; j < i; j++)
2911 if (scale_tot[i][j] != 0)
2917 hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2920 for (cg = cg0; cg < cg1; cg++)
2922 svmul(scale, cg_cm[cg], cgsc);
2923 a0 = cgs->index[cg];
2924 a1 = cgs->index[cg+1];
2925 for (a = a0; a < a1; a++)
2927 if (distance2(cgsc, x[a]) > hbuf2)
2936 /* Anistropic scaling */
2937 for (cg = cg0; cg < cg1; cg++)
2939 /* Since scale_tot contains the transpose of the scaling matrix,
2940 * we need to multiply with the transpose.
2942 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2943 a0 = cgs->index[cg];
2944 a1 = cgs->index[cg+1];
2945 for (a = a0; a < a1; a++)
2947 if (distance2(cgsc, x[a]) > hbuf2)