2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/utility/smalloc.h"
52 #include "mtop_util.h"
53 #include "gmx_omp_nthreads.h"
55 #include "gromacs/utility/gmxomp.h"
57 /* Routines to send/recieve coordinates and force
58 * of constructing atoms.
61 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
65 return pbc_dx_aiuc(pbc, xi, xj, dx);
74 /* Vsite construction routines */
76 static void constr_vsite2(rvec xi, rvec xj, rvec x, real a, t_pbc *pbc)
86 pbc_dx_aiuc(pbc, xj, xi, dx);
87 x[XX] = xi[XX] + a*dx[XX];
88 x[YY] = xi[YY] + a*dx[YY];
89 x[ZZ] = xi[ZZ] + a*dx[ZZ];
93 x[XX] = b*xi[XX] + a*xj[XX];
94 x[YY] = b*xi[YY] + a*xj[YY];
95 x[ZZ] = b*xi[ZZ] + a*xj[ZZ];
102 static void constr_vsite3(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
113 pbc_dx_aiuc(pbc, xj, xi, dxj);
114 pbc_dx_aiuc(pbc, xk, xi, dxk);
115 x[XX] = xi[XX] + a*dxj[XX] + b*dxk[XX];
116 x[YY] = xi[YY] + a*dxj[YY] + b*dxk[YY];
117 x[ZZ] = xi[ZZ] + a*dxj[ZZ] + b*dxk[ZZ];
121 x[XX] = c*xi[XX] + a*xj[XX] + b*xk[XX];
122 x[YY] = c*xi[YY] + a*xj[YY] + b*xk[YY];
123 x[ZZ] = c*xi[ZZ] + a*xj[ZZ] + b*xk[ZZ];
127 /* TOTAL: 17 flops */
130 static void constr_vsite3FD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
136 pbc_rvec_sub(pbc, xj, xi, xij);
137 pbc_rvec_sub(pbc, xk, xj, xjk);
140 /* temp goes from i to a point on the line jk */
141 temp[XX] = xij[XX] + a*xjk[XX];
142 temp[YY] = xij[YY] + a*xjk[YY];
143 temp[ZZ] = xij[ZZ] + a*xjk[ZZ];
146 c = b*gmx_invsqrt(iprod(temp, temp));
149 x[XX] = xi[XX] + c*temp[XX];
150 x[YY] = xi[YY] + c*temp[YY];
151 x[ZZ] = xi[ZZ] + c*temp[ZZ];
154 /* TOTAL: 34 flops */
157 static void constr_vsite3FAD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b, t_pbc *pbc)
160 real a1, b1, c1, invdij;
162 pbc_rvec_sub(pbc, xj, xi, xij);
163 pbc_rvec_sub(pbc, xk, xj, xjk);
166 invdij = gmx_invsqrt(iprod(xij, xij));
167 c1 = invdij * invdij * iprod(xij, xjk);
168 xp[XX] = xjk[XX] - c1*xij[XX];
169 xp[YY] = xjk[YY] - c1*xij[YY];
170 xp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
172 b1 = b*gmx_invsqrt(iprod(xp, xp));
175 x[XX] = xi[XX] + a1*xij[XX] + b1*xp[XX];
176 x[YY] = xi[YY] + a1*xij[YY] + b1*xp[YY];
177 x[ZZ] = xi[ZZ] + a1*xij[ZZ] + b1*xp[ZZ];
180 /* TOTAL: 63 flops */
183 static void constr_vsite3OUT(rvec xi, rvec xj, rvec xk, rvec x,
184 real a, real b, real c, t_pbc *pbc)
188 pbc_rvec_sub(pbc, xj, xi, xij);
189 pbc_rvec_sub(pbc, xk, xi, xik);
190 cprod(xij, xik, temp);
193 x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX];
194 x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY];
195 x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ];
198 /* TOTAL: 33 flops */
201 static void constr_vsite4FD(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
202 real a, real b, real c, t_pbc *pbc)
204 rvec xij, xjk, xjl, temp;
207 pbc_rvec_sub(pbc, xj, xi, xij);
208 pbc_rvec_sub(pbc, xk, xj, xjk);
209 pbc_rvec_sub(pbc, xl, xj, xjl);
212 /* temp goes from i to a point on the plane jkl */
213 temp[XX] = xij[XX] + a*xjk[XX] + b*xjl[XX];
214 temp[YY] = xij[YY] + a*xjk[YY] + b*xjl[YY];
215 temp[ZZ] = xij[ZZ] + a*xjk[ZZ] + b*xjl[ZZ];
218 d = c*gmx_invsqrt(iprod(temp, temp));
221 x[XX] = xi[XX] + d*temp[XX];
222 x[YY] = xi[YY] + d*temp[YY];
223 x[ZZ] = xi[ZZ] + d*temp[ZZ];
226 /* TOTAL: 43 flops */
230 static void constr_vsite4FDN(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
231 real a, real b, real c, t_pbc *pbc)
233 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
236 pbc_rvec_sub(pbc, xj, xi, xij);
237 pbc_rvec_sub(pbc, xk, xi, xik);
238 pbc_rvec_sub(pbc, xl, xi, xil);
251 rvec_sub(ra, xij, rja);
252 rvec_sub(rb, xij, rjb);
258 d = c*gmx_invsqrt(norm2(rm));
261 x[XX] = xi[XX] + d*rm[XX];
262 x[YY] = xi[YY] + d*rm[YY];
263 x[ZZ] = xi[ZZ] + d*rm[ZZ];
266 /* TOTAL: 47 flops */
270 static int constr_vsiten(t_iatom *ia, t_iparams ip[],
278 n3 = 3*ip[ia[0]].vsiten.n;
281 copy_rvec(x[ai], x1);
283 for (i = 3; i < n3; i += 3)
286 a = ip[ia[i]].vsiten.a;
289 pbc_dx_aiuc(pbc, x[ai], x1, dx);
293 rvec_sub(x[ai], x1, dx);
295 dsum[XX] += a*dx[XX];
296 dsum[YY] += a*dx[YY];
297 dsum[ZZ] += a*dx[ZZ];
301 x[av][XX] = x1[XX] + dsum[XX];
302 x[av][YY] = x1[YY] + dsum[YY];
303 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
309 void construct_vsites_thread(gmx_vsite_t *vsite,
312 t_iparams ip[], t_ilist ilist[],
316 rvec xpbc, xv, vv, dx;
317 real a1, b1, c1, inv_dt;
318 int i, inc, ii, nra, nr, tp, ftype;
319 t_iatom avsite, ai, aj, ak, al, pbc_atom;
322 int *vsite_pbc, ishift;
323 rvec reftmp, vtmp, rtmp;
334 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
338 for (ftype = 0; (ftype < F_NRE); ftype++)
340 if ((interaction_function[ftype].flags & IF_VSITE) &&
343 nra = interaction_function[ftype].nratoms;
345 nr = ilist[ftype].nr;
346 ia = ilist[ftype].iatoms;
350 pbc_null2 = pbc_null;
352 else if (pbc_null != NULL)
354 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
357 for (i = 0; i < nr; )
361 /* The vsite and constructing atoms */
365 /* Constants for constructing vsites */
367 /* Check what kind of pbc we need to use */
370 /* No charge groups, vsite follows its own pbc */
372 copy_rvec(x[avsite], xpbc);
374 else if (vsite_pbc != NULL)
376 pbc_atom = vsite_pbc[i/(1+nra)];
381 /* We need to copy the coordinates here,
382 * single for single atom cg's pbc_atom
383 * is the vsite itself.
385 copy_rvec(x[pbc_atom], xpbc);
387 pbc_null2 = pbc_null;
398 /* Copy the old position */
399 copy_rvec(x[avsite], xv);
401 /* Construct the vsite depending on type */
406 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
412 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
418 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
424 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
431 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
439 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
448 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
452 inc = constr_vsiten(ia, ip, x, pbc_null2);
455 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
456 ftype, __FILE__, __LINE__);
461 /* Match the pbc of this vsite to the rest of its charge group */
462 ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
463 if (ishift != CENTRAL)
465 rvec_add(xpbc, dx, x[avsite]);
470 /* Calculate velocity of vsite... */
471 rvec_sub(x[avsite], xv, vv);
472 svmul(inv_dt, vv, v[avsite]);
475 /* Increment loop variables */
483 void construct_vsites(gmx_vsite_t *vsite,
486 t_iparams ip[], t_ilist ilist[],
487 int ePBC, gmx_bool bMolPBC,
488 t_commrec *cr, matrix box)
490 t_pbc pbc, *pbc_null;
494 bDomDec = cr && DOMAINDECOMP(cr);
496 /* We only need to do pbc when we have inter-cg vsites */
497 if (ePBC != epbcNONE && (bDomDec || bMolPBC) && vsite->n_intercg_vsite)
499 /* This is wasting some CPU time as we now do this multiple times
500 * per MD step. But how often do we have vsites with full pbc?
502 pbc_null = set_pbc_dd(&pbc, ePBC, cr != NULL ? cr->dd : NULL, FALSE, box);
511 dd_move_x_vsites(cr->dd, box, x);
514 if (vsite->nthreads == 1)
516 construct_vsites_thread(vsite,
523 #pragma omp parallel num_threads(vsite->nthreads)
525 construct_vsites_thread(vsite,
527 ip, vsite->tdata[gmx_omp_get_thread_num()].ilist,
530 /* Now we can construct the vsites that might depend on other vsites */
531 construct_vsites_thread(vsite,
533 ip, vsite->tdata[vsite->nthreads].ilist,
538 static void spread_vsite2(t_iatom ia[], real a,
539 rvec x[], rvec f[], rvec fshift[],
540 t_pbc *pbc, t_graph *g)
552 svmul(1-a, f[av], fi);
553 svmul( a, f[av], fj);
562 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
564 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
569 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
570 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
578 if (fshift && (siv != CENTRAL || sij != CENTRAL))
580 rvec_inc(fshift[siv], f[av]);
581 rvec_dec(fshift[CENTRAL], fi);
582 rvec_dec(fshift[sij], fj);
585 /* TOTAL: 13 flops */
588 void construct_vsites_mtop(gmx_vsite_t *vsite,
589 gmx_mtop_t *mtop, rvec x[])
592 gmx_molblock_t *molb;
596 for (mb = 0; mb < mtop->nmolblock; mb++)
598 molb = &mtop->molblock[mb];
599 molt = &mtop->moltype[molb->type];
600 for (mol = 0; mol < molb->nmol; mol++)
602 construct_vsites(vsite, x+as, 0.0, NULL,
603 mtop->ffparams.iparams, molt->ilist,
604 epbcNONE, TRUE, NULL, NULL);
605 as += molt->atoms.nr;
610 static void spread_vsite3(t_iatom ia[], real a, real b,
611 rvec x[], rvec f[], rvec fshift[],
612 t_pbc *pbc, t_graph *g)
615 atom_id av, ai, aj, ak;
624 svmul(1-a-b, f[av], fi);
625 svmul( a, f[av], fj);
626 svmul( b, f[av], fk);
636 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
638 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
640 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
645 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
646 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
647 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
656 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
658 rvec_inc(fshift[siv], f[av]);
659 rvec_dec(fshift[CENTRAL], fi);
660 rvec_dec(fshift[sij], fj);
661 rvec_dec(fshift[sik], fk);
664 /* TOTAL: 20 flops */
667 static void spread_vsite3FD(t_iatom ia[], real a, real b,
668 rvec x[], rvec f[], rvec fshift[],
669 gmx_bool VirCorr, matrix dxdf,
670 t_pbc *pbc, t_graph *g)
672 real fx, fy, fz, c, invl, fproj, a1;
673 rvec xvi, xij, xjk, xix, fv, temp;
674 t_iatom av, ai, aj, ak;
675 int svi, sji, skj, d;
682 copy_rvec(f[av], fv);
684 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
685 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
688 /* xix goes from i to point x on the line jk */
689 xix[XX] = xij[XX]+a*xjk[XX];
690 xix[YY] = xij[YY]+a*xjk[YY];
691 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
694 invl = gmx_invsqrt(iprod(xix, xix));
698 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
700 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
701 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
702 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
705 /* c is already calculated in constr_vsite3FD
706 storing c somewhere will save 26 flops! */
709 f[ai][XX] += fv[XX] - temp[XX];
710 f[ai][YY] += fv[YY] - temp[YY];
711 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
712 f[aj][XX] += a1*temp[XX];
713 f[aj][YY] += a1*temp[YY];
714 f[aj][ZZ] += a1*temp[ZZ];
715 f[ak][XX] += a*temp[XX];
716 f[ak][YY] += a*temp[YY];
717 f[ak][ZZ] += a*temp[ZZ];
722 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
724 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
726 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
731 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
738 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
740 rvec_dec(fshift[svi], fv);
741 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
742 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
743 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
744 fshift[ sji][XX] += temp[XX];
745 fshift[ sji][YY] += temp[YY];
746 fshift[ sji][ZZ] += temp[ZZ];
747 fshift[ skj][XX] += a*temp[XX];
748 fshift[ skj][YY] += a*temp[YY];
749 fshift[ skj][ZZ] += a*temp[ZZ];
754 /* When VirCorr=TRUE, the virial for the current forces is not
755 * calculated from the redistributed forces. This means that
756 * the effect of non-linear virtual site constructions on the virial
757 * needs to be added separately. This contribution can be calculated
758 * in many ways, but the simplest and cheapest way is to use
759 * the first constructing atom ai as a reference position in space:
760 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
765 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
767 for (i = 0; i < DIM; i++)
769 for (j = 0; j < DIM; j++)
771 /* As xix is a linear combination of j and k, use that here */
772 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
777 /* TOTAL: 61 flops */
780 static void spread_vsite3FAD(t_iatom ia[], real a, real b,
781 rvec x[], rvec f[], rvec fshift[],
782 gmx_bool VirCorr, matrix dxdf,
783 t_pbc *pbc, t_graph *g)
785 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
786 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
787 t_iatom av, ai, aj, ak;
788 int svi, sji, skj, d;
795 copy_rvec(f[ia[1]], fv);
797 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
798 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
801 invdij = gmx_invsqrt(iprod(xij, xij));
802 invdij2 = invdij * invdij;
803 c1 = iprod(xij, xjk) * invdij2;
804 xperp[XX] = xjk[XX] - c1*xij[XX];
805 xperp[YY] = xjk[YY] - c1*xij[YY];
806 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
807 /* xperp in plane ijk, perp. to ij */
808 invdp = gmx_invsqrt(iprod(xperp, xperp));
813 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
814 storing them somewhere will save 45 flops! */
816 fproj = iprod(xij, fv)*invdij2;
817 svmul(fproj, xij, Fpij); /* proj. f on xij */
818 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
819 svmul(b1*fproj, xperp, f3);
822 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
823 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
824 for (d = 0; (d < DIM); d++)
832 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
833 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
834 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
835 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
836 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
837 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
845 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
847 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
849 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
854 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
861 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
863 rvec_dec(fshift[svi], fv);
864 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
865 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
866 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
867 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
868 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
869 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
870 fshift[ skj][XX] += f2[XX];
871 fshift[ skj][YY] += f2[YY];
872 fshift[ skj][ZZ] += f2[ZZ];
880 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
882 for (i = 0; i < DIM; i++)
884 for (j = 0; j < DIM; j++)
886 /* Note that xik=xij+xjk, so we have to add xij*f2 */
889 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
895 /* TOTAL: 113 flops */
898 static void spread_vsite3OUT(t_iatom ia[], real a, real b, real c,
899 rvec x[], rvec f[], rvec fshift[],
900 gmx_bool VirCorr, matrix dxdf,
901 t_pbc *pbc, t_graph *g)
903 rvec xvi, xij, xik, fv, fj, fk;
905 atom_id av, ai, aj, ak;
914 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
915 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
918 copy_rvec(f[av], fv);
925 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
926 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
927 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
929 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
930 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
931 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
934 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
935 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
936 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
943 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
945 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
947 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
952 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
959 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
961 rvec_dec(fshift[svi], fv);
962 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
963 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
964 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
965 rvec_inc(fshift[sji], fj);
966 rvec_inc(fshift[ski], fk);
974 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
976 for (i = 0; i < DIM; i++)
978 for (j = 0; j < DIM; j++)
980 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
985 /* TOTAL: 54 flops */
988 static void spread_vsite4FD(t_iatom ia[], real a, real b, real c,
989 rvec x[], rvec f[], rvec fshift[],
990 gmx_bool VirCorr, matrix dxdf,
991 t_pbc *pbc, t_graph *g)
993 real d, invl, fproj, a1;
994 rvec xvi, xij, xjk, xjl, xix, fv, temp;
995 atom_id av, ai, aj, ak, al;
997 int svi, sji, skj, slj, m;
1005 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1006 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1007 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1010 /* xix goes from i to point x on the plane jkl */
1011 for (m = 0; m < DIM; m++)
1013 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1017 invl = gmx_invsqrt(iprod(xix, xix));
1019 /* 4 + ?10? flops */
1021 copy_rvec(f[av], fv);
1023 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1025 for (m = 0; m < DIM; m++)
1027 temp[m] = d*(fv[m] - fproj*xix[m]);
1031 /* c is already calculated in constr_vsite3FD
1032 storing c somewhere will save 35 flops! */
1035 for (m = 0; m < DIM; m++)
1037 f[ai][m] += fv[m] - temp[m];
1038 f[aj][m] += a1*temp[m];
1039 f[ak][m] += a*temp[m];
1040 f[al][m] += b*temp[m];
1046 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1048 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1050 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1052 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1057 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1065 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1067 rvec_dec(fshift[svi], fv);
1068 for (m = 0; m < DIM; m++)
1070 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1071 fshift[ sji][m] += temp[m];
1072 fshift[ skj][m] += a*temp[m];
1073 fshift[ slj][m] += b*temp[m];
1082 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1084 for (i = 0; i < DIM; i++)
1086 for (j = 0; j < DIM; j++)
1088 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1093 /* TOTAL: 77 flops */
1097 static void spread_vsite4FDN(t_iatom ia[], real a, real b, real c,
1098 rvec x[], rvec f[], rvec fshift[],
1099 gmx_bool VirCorr, matrix dxdf,
1100 t_pbc *pbc, t_graph *g)
1102 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1103 rvec fv, fj, fk, fl;
1107 int av, ai, aj, ak, al;
1108 int svi, sij, sik, sil;
1110 /* DEBUG: check atom indices */
1117 copy_rvec(f[av], fv);
1119 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1120 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1121 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1134 rvec_sub(ra, xij, rja);
1135 rvec_sub(rb, xij, rjb);
1136 rvec_sub(rb, ra, rab);
1139 cprod(rja, rjb, rm);
1142 invrm = gmx_invsqrt(norm2(rm));
1143 denom = invrm*invrm;
1146 cfx = c*invrm*fv[XX];
1147 cfy = c*invrm*fv[YY];
1148 cfz = c*invrm*fv[ZZ];
1159 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1160 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1161 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1172 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1173 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1174 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1185 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1186 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1187 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1190 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1191 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1192 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1193 rvec_inc(f[aj], fj);
1194 rvec_inc(f[ak], fk);
1195 rvec_inc(f[al], fl);
1200 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1202 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1204 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1206 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1211 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1218 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1220 rvec_dec(fshift[svi], fv);
1221 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1222 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1223 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1224 rvec_inc(fshift[sij], fj);
1225 rvec_inc(fshift[sik], fk);
1226 rvec_inc(fshift[sil], fl);
1234 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1236 for (i = 0; i < DIM; i++)
1238 for (j = 0; j < DIM; j++)
1240 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1245 /* Total: 207 flops (Yuck!) */
1249 static int spread_vsiten(t_iatom ia[], t_iparams ip[],
1250 rvec x[], rvec f[], rvec fshift[],
1251 t_pbc *pbc, t_graph *g)
1259 n3 = 3*ip[ia[0]].vsiten.n;
1261 copy_rvec(x[av], xv);
1263 for (i = 0; i < n3; i += 3)
1268 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1273 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1279 a = ip[ia[i]].vsiten.a;
1280 svmul(a, f[av], fi);
1281 rvec_inc(f[ai], fi);
1282 if (fshift && siv != CENTRAL)
1284 rvec_inc(fshift[siv], fi);
1285 rvec_dec(fshift[CENTRAL], fi);
1294 static int vsite_count(const t_ilist *ilist, int ftype)
1296 if (ftype == F_VSITEN)
1298 return ilist[ftype].nr/3;
1302 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1306 static void spread_vsite_f_thread(gmx_vsite_t *vsite,
1307 rvec x[], rvec f[], rvec *fshift,
1308 gmx_bool VirCorr, matrix dxdf,
1309 t_iparams ip[], t_ilist ilist[],
1310 t_graph *g, t_pbc *pbc_null)
1314 int i, inc, m, nra, nr, tp, ftype;
1324 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
1326 /* this loop goes backwards to be able to build *
1327 * higher type vsites from lower types */
1330 for (ftype = F_NRE-1; (ftype >= 0); ftype--)
1332 if ((interaction_function[ftype].flags & IF_VSITE) &&
1333 ilist[ftype].nr > 0)
1335 nra = interaction_function[ftype].nratoms;
1337 nr = ilist[ftype].nr;
1338 ia = ilist[ftype].iatoms;
1342 pbc_null2 = pbc_null;
1344 else if (pbc_null != NULL)
1346 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
1349 for (i = 0; i < nr; )
1351 if (vsite_pbc != NULL)
1353 if (vsite_pbc[i/(1+nra)] > -2)
1355 pbc_null2 = pbc_null;
1365 /* Constants for constructing */
1366 a1 = ip[tp].vsite.a;
1367 /* Construct the vsite depending on type */
1371 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1374 b1 = ip[tp].vsite.b;
1375 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1378 b1 = ip[tp].vsite.b;
1379 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1382 b1 = ip[tp].vsite.b;
1383 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1386 b1 = ip[tp].vsite.b;
1387 c1 = ip[tp].vsite.c;
1388 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1391 b1 = ip[tp].vsite.b;
1392 c1 = ip[tp].vsite.c;
1393 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1396 b1 = ip[tp].vsite.b;
1397 c1 = ip[tp].vsite.c;
1398 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1401 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1404 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1405 ftype, __FILE__, __LINE__);
1407 clear_rvec(f[ia[1]]);
1409 /* Increment loop variables */
1417 void spread_vsite_f(gmx_vsite_t *vsite,
1418 rvec x[], rvec f[], rvec *fshift,
1419 gmx_bool VirCorr, matrix vir,
1420 t_nrnb *nrnb, t_idef *idef,
1421 int ePBC, gmx_bool bMolPBC, t_graph *g, matrix box,
1424 t_pbc pbc, *pbc_null;
1427 /* We only need to do pbc when we have inter-cg vsites */
1428 if ((DOMAINDECOMP(cr) || bMolPBC) && vsite->n_intercg_vsite)
1430 /* This is wasting some CPU time as we now do this multiple times
1431 * per MD step. But how often do we have vsites with full pbc?
1433 pbc_null = set_pbc_dd(&pbc, ePBC, cr->dd, FALSE, box);
1440 if (DOMAINDECOMP(cr))
1442 dd_clear_f_vsites(cr->dd, f);
1445 if (vsite->nthreads == 1)
1447 spread_vsite_f_thread(vsite,
1449 VirCorr, vsite->tdata[0].dxdf,
1450 idef->iparams, idef->il,
1455 /* First spread the vsites that might depend on other vsites */
1456 spread_vsite_f_thread(vsite,
1458 VirCorr, vsite->tdata[vsite->nthreads].dxdf,
1460 vsite->tdata[vsite->nthreads].ilist,
1463 #pragma omp parallel num_threads(vsite->nthreads)
1468 thread = gmx_omp_get_thread_num();
1470 if (thread == 0 || fshift == NULL)
1478 fshift_t = vsite->tdata[thread].fshift;
1480 for (i = 0; i < SHIFTS; i++)
1482 clear_rvec(fshift_t[i]);
1486 spread_vsite_f_thread(vsite,
1488 VirCorr, vsite->tdata[thread].dxdf,
1490 vsite->tdata[thread].ilist,
1498 for (th = 1; th < vsite->nthreads; th++)
1500 for (i = 0; i < SHIFTS; i++)
1502 rvec_inc(fshift[i], vsite->tdata[th].fshift[i]);
1512 for (th = 0; th < (vsite->nthreads == 1 ? 1 : vsite->nthreads+1); th++)
1514 for (i = 0; i < DIM; i++)
1516 for (j = 0; j < DIM; j++)
1518 vir[i][j] += -0.5*vsite->tdata[th].dxdf[i][j];
1524 if (DOMAINDECOMP(cr))
1526 dd_move_f_vsites(cr->dd, f, fshift);
1529 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1530 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1531 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1532 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1533 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1534 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1535 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1536 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1539 static int *atom2cg(t_block *cgs)
1543 snew(a2cg, cgs->index[cgs->nr]);
1544 for (cg = 0; cg < cgs->nr; cg++)
1546 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
1555 static int count_intercg_vsite(gmx_mtop_t *mtop,
1556 gmx_bool *bHaveChargeGroups)
1558 int mb, mt, ftype, nral, i, cg, a;
1559 gmx_molblock_t *molb;
1560 gmx_moltype_t *molt;
1564 int n_intercg_vsite;
1566 *bHaveChargeGroups = FALSE;
1568 n_intercg_vsite = 0;
1569 for (mb = 0; mb < mtop->nmolblock; mb++)
1571 molb = &mtop->molblock[mb];
1572 molt = &mtop->moltype[molb->type];
1574 if (molt->cgs.nr < molt->atoms.nr)
1576 *bHaveChargeGroups = TRUE;
1579 a2cg = atom2cg(&molt->cgs);
1580 for (ftype = 0; ftype < F_NRE; ftype++)
1582 if (interaction_function[ftype].flags & IF_VSITE)
1585 il = &molt->ilist[ftype];
1587 for (i = 0; i < il->nr; i += 1+nral)
1590 for (a = 1; a < nral; a++)
1592 if (a2cg[ia[1+a]] != cg)
1594 n_intercg_vsite += molb->nmol;
1604 return n_intercg_vsite;
1607 static int **get_vsite_pbc(t_iparams *iparams, t_ilist *ilist,
1608 t_atom *atom, t_mdatoms *md,
1609 t_block *cgs, int *a2cg)
1611 int ftype, nral, i, j, vsi, vsite, cg_v, cg_c, a, nc3 = 0;
1614 int **vsite_pbc, *vsite_pbc_f;
1616 gmx_bool bViteOnlyCG_and_FirstAtom;
1618 /* Make an array that tells if the pbc of an atom is set */
1619 snew(pbc_set, cgs->index[cgs->nr]);
1620 /* PBC is set for all non vsites */
1621 for (a = 0; a < cgs->index[cgs->nr]; a++)
1623 if ((atom && atom[a].ptype != eptVSite) ||
1624 (md && md->ptype[a] != eptVSite))
1630 snew(vsite_pbc, F_VSITEN-F_VSITE2+1);
1632 for (ftype = 0; ftype < F_NRE; ftype++)
1634 if (interaction_function[ftype].flags & IF_VSITE)
1640 snew(vsite_pbc[ftype-F_VSITE2], il->nr/(1+nral));
1641 vsite_pbc_f = vsite_pbc[ftype-F_VSITE2];
1649 /* A value of -2 signals that this vsite and its contructing
1650 * atoms are all within the same cg, so no pbc is required.
1652 vsite_pbc_f[vsi] = -2;
1653 /* Check if constructing atoms are outside the vsite's cg */
1655 if (ftype == F_VSITEN)
1657 nc3 = 3*iparams[ia[i]].vsiten.n;
1658 for (j = 0; j < nc3; j += 3)
1660 if (a2cg[ia[i+j+2]] != cg_v)
1662 vsite_pbc_f[vsi] = -1;
1668 for (a = 1; a < nral; a++)
1670 if (a2cg[ia[i+1+a]] != cg_v)
1672 vsite_pbc_f[vsi] = -1;
1676 if (vsite_pbc_f[vsi] == -1)
1678 /* Check if this is the first processed atom of a vsite only cg */
1679 bViteOnlyCG_and_FirstAtom = TRUE;
1680 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1682 /* Non-vsites already have pbc set, so simply check for pbc_set */
1685 bViteOnlyCG_and_FirstAtom = FALSE;
1689 if (bViteOnlyCG_and_FirstAtom)
1691 /* First processed atom of a vsite only charge group.
1692 * The pbc of the input coordinates to construct_vsites
1693 * should be preserved.
1695 vsite_pbc_f[vsi] = vsite;
1697 else if (cg_v != a2cg[ia[1+i+1]])
1699 /* This vsite has a different charge group index
1700 * than it's first constructing atom
1701 * and the charge group has more than one atom,
1702 * search for the first normal particle
1703 * or vsite that already had its pbc defined.
1704 * If nothing is found, use full pbc for this vsite.
1706 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1708 if (a != vsite && pbc_set[a])
1710 vsite_pbc_f[vsi] = a;
1713 fprintf(debug, "vsite %d match pbc with atom %d\n",
1721 fprintf(debug, "vsite atom %d cg %d - %d pbc atom %d\n",
1722 vsite+1, cgs->index[cg_v]+1, cgs->index[cg_v+1],
1723 vsite_pbc_f[vsi]+1);
1727 if (ftype == F_VSITEN)
1729 /* The other entries in vsite_pbc_f are not used for center vsites */
1737 /* This vsite now has its pbc defined */
1749 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop, t_commrec *cr,
1750 gmx_bool bSerial_NoPBC)
1756 gmx_moltype_t *molt;
1759 /* check if there are vsites */
1761 for (i = 0; i < F_NRE; i++)
1763 if (interaction_function[i].flags & IF_VSITE)
1765 nvsite += gmx_mtop_ftype_count(mtop, i);
1776 vsite->n_intercg_vsite = count_intercg_vsite(mtop,
1777 &vsite->bHaveChargeGroups);
1779 /* If we don't have charge groups, the vsite follows its own pbc */
1780 if (!bSerial_NoPBC &&
1781 vsite->bHaveChargeGroups &&
1782 vsite->n_intercg_vsite > 0 && DOMAINDECOMP(cr))
1784 vsite->nvsite_pbc_molt = mtop->nmoltype;
1785 snew(vsite->vsite_pbc_molt, vsite->nvsite_pbc_molt);
1786 for (mt = 0; mt < mtop->nmoltype; mt++)
1788 molt = &mtop->moltype[mt];
1789 /* Make an atom to charge group index */
1790 a2cg = atom2cg(&molt->cgs);
1791 vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop->ffparams.iparams,
1793 molt->atoms.atom, NULL,
1798 snew(vsite->vsite_pbc_loc_nalloc, F_VSITEN-F_VSITE2+1);
1799 snew(vsite->vsite_pbc_loc, F_VSITEN-F_VSITE2+1);
1804 vsite->nthreads = 1;
1808 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
1812 /* We need one extra thread data structure for the overlap vsites */
1813 snew(vsite->tdata, vsite->nthreads+1);
1816 vsite->th_ind = NULL;
1817 vsite->th_ind_nalloc = 0;
1822 static void prepare_vsite_thread(const t_ilist *ilist,
1823 gmx_vsite_thread_t *vsite_th)
1827 for (ftype = 0; ftype < F_NRE; ftype++)
1829 if (interaction_function[ftype].flags & IF_VSITE)
1831 if (ilist[ftype].nr > vsite_th->ilist[ftype].nalloc)
1833 vsite_th->ilist[ftype].nalloc = over_alloc_large(ilist[ftype].nr);
1834 srenew(vsite_th->ilist[ftype].iatoms, vsite_th->ilist[ftype].nalloc);
1837 vsite_th->ilist[ftype].nr = 0;
1842 void split_vsites_over_threads(const t_ilist *ilist,
1843 const t_iparams *ip,
1844 const t_mdatoms *mdatoms,
1845 gmx_bool bLimitRange,
1849 int vsite_atom_range, natperthread;
1854 int nral1, inc, i, j;
1856 if (vsite->nthreads == 1)
1862 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
1863 for (th = 0; th < vsite->nthreads; th++)
1865 prepare_vsite_thread(ilist, &vsite->tdata[th]);
1867 /* Master threads does the (potential) overlap vsites */
1868 prepare_vsite_thread(ilist, &vsite->tdata[vsite->nthreads]);
1870 /* The current way of distributing the vsites over threads in primitive.
1871 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
1872 * without taking into account how the vsites are distributed.
1873 * Without domain decomposition we bLimitRange=TRUE and we at least
1874 * tighten the upper bound of the range (useful for common systems
1875 * such as a vsite-protein in 3-site water).
1879 vsite_atom_range = -1;
1880 for (ftype = 0; ftype < F_NRE; ftype++)
1882 if (interaction_function[ftype].flags & IF_VSITE)
1884 if (ftype != F_VSITEN)
1886 nral1 = 1 + NRAL(ftype);
1887 iat = ilist[ftype].iatoms;
1888 for (i = 0; i < ilist[ftype].nr; i += nral1)
1890 for (j = i + 1; j < i + nral1; j++)
1892 vsite_atom_range = max(vsite_atom_range, iat[j]);
1900 iat = ilist[ftype].iatoms;
1903 while (i < ilist[ftype].nr)
1905 /* The 3 below is from 1+NRAL(ftype)=3 */
1906 vs_ind_end = i + ip[iat[i]].vsiten.n*3;
1908 vsite_atom_range = max(vsite_atom_range, iat[i+1]);
1909 while (i < vs_ind_end)
1911 vsite_atom_range = max(vsite_atom_range, iat[i+2]);
1922 vsite_atom_range = mdatoms->homenr;
1924 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
1928 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
1931 /* To simplify the vsite assignment, we make an index which tells us
1932 * to which thread particles, both non-vsites and vsites, are assigned.
1934 if (mdatoms->nr > vsite->th_ind_nalloc)
1936 vsite->th_ind_nalloc = over_alloc_large(mdatoms->nr);
1937 srenew(vsite->th_ind, vsite->th_ind_nalloc);
1939 th_ind = vsite->th_ind;
1941 for (i = 0; i < mdatoms->nr; i++)
1943 if (mdatoms->ptype[i] == eptVSite)
1945 /* vsites are not assigned to a thread yet */
1950 /* assign non-vsite particles to thread th */
1953 if (i == (th + 1)*natperthread && th < vsite->nthreads)
1959 for (ftype = 0; ftype < F_NRE; ftype++)
1961 if (interaction_function[ftype].flags & IF_VSITE)
1963 nral1 = 1 + NRAL(ftype);
1965 iat = ilist[ftype].iatoms;
1966 for (i = 0; i < ilist[ftype].nr; )
1968 th = iat[1+i]/natperthread;
1969 /* We would like to assign this vsite the thread th,
1970 * but it might depend on atoms outside the atom range of th
1971 * or on another vsite not assigned to thread th.
1973 if (ftype != F_VSITEN)
1975 for (j = i + 2; j < i + nral1; j++)
1977 if (th_ind[iat[j]] != th)
1979 /* Some constructing atoms are not assigned to
1980 * thread th, move this vsite to a separate batch.
1982 th = vsite->nthreads;
1988 /* The 3 below is from 1+NRAL(ftype)=3 */
1989 inc = ip[iat[i]].vsiten.n*3;
1990 for (j = i + 2; j < i + inc; j += 3)
1992 if (th_ind[iat[j]] != th)
1994 th = vsite->nthreads;
1998 /* Copy this vsite to the thread data struct of thread th */
1999 il_th = &vsite->tdata[th].ilist[ftype];
2000 for (j = i; j < i + inc; j++)
2002 il_th->iatoms[il_th->nr++] = iat[j];
2004 /* Update this vsite's thread index entry */
2005 th_ind[iat[1+i]] = th;
2014 for (ftype = 0; ftype < F_NRE; ftype++)
2016 if ((interaction_function[ftype].flags & IF_VSITE) &&
2017 ilist[ftype].nr > 0)
2019 fprintf(debug, "%-20s thread dist:",
2020 interaction_function[ftype].longname);
2021 for (th = 0; th < vsite->nthreads+1; th++)
2023 fprintf(debug, " %4d", vsite->tdata[th].ilist[ftype].nr);
2025 fprintf(debug, "\n");
2031 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
2036 if (vsite->n_intercg_vsite > 0)
2038 if (vsite->bHaveChargeGroups)
2040 /* Make an atom to charge group index */
2041 a2cg = atom2cg(&top->cgs);
2042 vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2043 top->idef.il, NULL, md,
2050 if (vsite->nthreads > 1)
2052 if (vsite->bHaveChargeGroups)
2054 gmx_fatal(FARGS, "The combination of threading, virtual sites and charge groups is not implemented");
2057 split_vsites_over_threads(top->idef.il, top->idef.iparams,
2058 md, !DOMAINDECOMP(cr), vsite);