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.
44 #include "types/commrec.h"
47 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/math/vec.h"
54 #include "mtop_util.h"
55 #include "gmx_omp_nthreads.h"
57 #include "gromacs/utility/gmxomp.h"
59 /* Routines to send/recieve coordinates and force
60 * of constructing atoms.
63 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
67 return pbc_dx_aiuc(pbc, xi, xj, dx);
76 /* Vsite construction routines */
78 static void constr_vsite2(rvec xi, rvec xj, rvec x, real a, t_pbc *pbc)
88 pbc_dx_aiuc(pbc, xj, xi, dx);
89 x[XX] = xi[XX] + a*dx[XX];
90 x[YY] = xi[YY] + a*dx[YY];
91 x[ZZ] = xi[ZZ] + a*dx[ZZ];
95 x[XX] = b*xi[XX] + a*xj[XX];
96 x[YY] = b*xi[YY] + a*xj[YY];
97 x[ZZ] = b*xi[ZZ] + a*xj[ZZ];
101 /* TOTAL: 10 flops */
104 static void constr_vsite3(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
115 pbc_dx_aiuc(pbc, xj, xi, dxj);
116 pbc_dx_aiuc(pbc, xk, xi, dxk);
117 x[XX] = xi[XX] + a*dxj[XX] + b*dxk[XX];
118 x[YY] = xi[YY] + a*dxj[YY] + b*dxk[YY];
119 x[ZZ] = xi[ZZ] + a*dxj[ZZ] + b*dxk[ZZ];
123 x[XX] = c*xi[XX] + a*xj[XX] + b*xk[XX];
124 x[YY] = c*xi[YY] + a*xj[YY] + b*xk[YY];
125 x[ZZ] = c*xi[ZZ] + a*xj[ZZ] + b*xk[ZZ];
129 /* TOTAL: 17 flops */
132 static void constr_vsite3FD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b,
138 pbc_rvec_sub(pbc, xj, xi, xij);
139 pbc_rvec_sub(pbc, xk, xj, xjk);
142 /* temp goes from i to a point on the line jk */
143 temp[XX] = xij[XX] + a*xjk[XX];
144 temp[YY] = xij[YY] + a*xjk[YY];
145 temp[ZZ] = xij[ZZ] + a*xjk[ZZ];
148 c = b*gmx_invsqrt(iprod(temp, temp));
151 x[XX] = xi[XX] + c*temp[XX];
152 x[YY] = xi[YY] + c*temp[YY];
153 x[ZZ] = xi[ZZ] + c*temp[ZZ];
156 /* TOTAL: 34 flops */
159 static void constr_vsite3FAD(rvec xi, rvec xj, rvec xk, rvec x, real a, real b, t_pbc *pbc)
162 real a1, b1, c1, invdij;
164 pbc_rvec_sub(pbc, xj, xi, xij);
165 pbc_rvec_sub(pbc, xk, xj, xjk);
168 invdij = gmx_invsqrt(iprod(xij, xij));
169 c1 = invdij * invdij * iprod(xij, xjk);
170 xp[XX] = xjk[XX] - c1*xij[XX];
171 xp[YY] = xjk[YY] - c1*xij[YY];
172 xp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
174 b1 = b*gmx_invsqrt(iprod(xp, xp));
177 x[XX] = xi[XX] + a1*xij[XX] + b1*xp[XX];
178 x[YY] = xi[YY] + a1*xij[YY] + b1*xp[YY];
179 x[ZZ] = xi[ZZ] + a1*xij[ZZ] + b1*xp[ZZ];
182 /* TOTAL: 63 flops */
185 static void constr_vsite3OUT(rvec xi, rvec xj, rvec xk, rvec x,
186 real a, real b, real c, t_pbc *pbc)
190 pbc_rvec_sub(pbc, xj, xi, xij);
191 pbc_rvec_sub(pbc, xk, xi, xik);
192 cprod(xij, xik, temp);
195 x[XX] = xi[XX] + a*xij[XX] + b*xik[XX] + c*temp[XX];
196 x[YY] = xi[YY] + a*xij[YY] + b*xik[YY] + c*temp[YY];
197 x[ZZ] = xi[ZZ] + a*xij[ZZ] + b*xik[ZZ] + c*temp[ZZ];
200 /* TOTAL: 33 flops */
203 static void constr_vsite4FD(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
204 real a, real b, real c, t_pbc *pbc)
206 rvec xij, xjk, xjl, temp;
209 pbc_rvec_sub(pbc, xj, xi, xij);
210 pbc_rvec_sub(pbc, xk, xj, xjk);
211 pbc_rvec_sub(pbc, xl, xj, xjl);
214 /* temp goes from i to a point on the plane jkl */
215 temp[XX] = xij[XX] + a*xjk[XX] + b*xjl[XX];
216 temp[YY] = xij[YY] + a*xjk[YY] + b*xjl[YY];
217 temp[ZZ] = xij[ZZ] + a*xjk[ZZ] + b*xjl[ZZ];
220 d = c*gmx_invsqrt(iprod(temp, temp));
223 x[XX] = xi[XX] + d*temp[XX];
224 x[YY] = xi[YY] + d*temp[YY];
225 x[ZZ] = xi[ZZ] + d*temp[ZZ];
228 /* TOTAL: 43 flops */
232 static void constr_vsite4FDN(rvec xi, rvec xj, rvec xk, rvec xl, rvec x,
233 real a, real b, real c, t_pbc *pbc)
235 rvec xij, xik, xil, ra, rb, rja, rjb, rm;
238 pbc_rvec_sub(pbc, xj, xi, xij);
239 pbc_rvec_sub(pbc, xk, xi, xik);
240 pbc_rvec_sub(pbc, xl, xi, xil);
253 rvec_sub(ra, xij, rja);
254 rvec_sub(rb, xij, rjb);
260 d = c*gmx_invsqrt(norm2(rm));
263 x[XX] = xi[XX] + d*rm[XX];
264 x[YY] = xi[YY] + d*rm[YY];
265 x[ZZ] = xi[ZZ] + d*rm[ZZ];
268 /* TOTAL: 47 flops */
272 static int constr_vsiten(t_iatom *ia, t_iparams ip[],
280 n3 = 3*ip[ia[0]].vsiten.n;
283 copy_rvec(x[ai], x1);
285 for (i = 3; i < n3; i += 3)
288 a = ip[ia[i]].vsiten.a;
291 pbc_dx_aiuc(pbc, x[ai], x1, dx);
295 rvec_sub(x[ai], x1, dx);
297 dsum[XX] += a*dx[XX];
298 dsum[YY] += a*dx[YY];
299 dsum[ZZ] += a*dx[ZZ];
303 x[av][XX] = x1[XX] + dsum[XX];
304 x[av][YY] = x1[YY] + dsum[YY];
305 x[av][ZZ] = x1[ZZ] + dsum[ZZ];
311 void construct_vsites_thread(gmx_vsite_t *vsite,
314 t_iparams ip[], t_ilist ilist[],
318 rvec xpbc, xv, vv, dx;
319 real a1, b1, c1, inv_dt;
320 int i, inc, ii, nra, nr, tp, ftype;
321 t_iatom avsite, ai, aj, ak, al, pbc_atom;
324 int *vsite_pbc, ishift;
325 rvec reftmp, vtmp, rtmp;
336 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
340 for (ftype = 0; (ftype < F_NRE); ftype++)
342 if ((interaction_function[ftype].flags & IF_VSITE) &&
345 nra = interaction_function[ftype].nratoms;
347 nr = ilist[ftype].nr;
348 ia = ilist[ftype].iatoms;
352 pbc_null2 = pbc_null;
354 else if (pbc_null != NULL)
356 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
359 for (i = 0; i < nr; )
363 /* The vsite and constructing atoms */
367 /* Constants for constructing vsites */
369 /* Check what kind of pbc we need to use */
372 /* No charge groups, vsite follows its own pbc */
374 copy_rvec(x[avsite], xpbc);
376 else if (vsite_pbc != NULL)
378 pbc_atom = vsite_pbc[i/(1+nra)];
383 /* We need to copy the coordinates here,
384 * single for single atom cg's pbc_atom
385 * is the vsite itself.
387 copy_rvec(x[pbc_atom], xpbc);
389 pbc_null2 = pbc_null;
400 /* Copy the old position */
401 copy_rvec(x[avsite], xv);
403 /* Construct the vsite depending on type */
408 constr_vsite2(x[ai], x[aj], x[avsite], a1, pbc_null2);
414 constr_vsite3(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
420 constr_vsite3FD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
426 constr_vsite3FAD(x[ai], x[aj], x[ak], x[avsite], a1, b1, pbc_null2);
433 constr_vsite3OUT(x[ai], x[aj], x[ak], x[avsite], a1, b1, c1, pbc_null2);
441 constr_vsite4FD(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
450 constr_vsite4FDN(x[ai], x[aj], x[ak], x[al], x[avsite], a1, b1, c1,
454 inc = constr_vsiten(ia, ip, x, pbc_null2);
457 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
458 ftype, __FILE__, __LINE__);
463 /* Match the pbc of this vsite to the rest of its charge group */
464 ishift = pbc_dx_aiuc(pbc_null, x[avsite], xpbc, dx);
465 if (ishift != CENTRAL)
467 rvec_add(xpbc, dx, x[avsite]);
472 /* Calculate velocity of vsite... */
473 rvec_sub(x[avsite], xv, vv);
474 svmul(inv_dt, vv, v[avsite]);
477 /* Increment loop variables */
485 void construct_vsites(gmx_vsite_t *vsite,
488 t_iparams ip[], t_ilist ilist[],
489 int ePBC, gmx_bool bMolPBC,
490 t_commrec *cr, matrix box)
492 t_pbc pbc, *pbc_null;
496 bDomDec = cr && DOMAINDECOMP(cr);
498 /* We only need to do pbc when we have inter-cg vsites */
499 if (ePBC != epbcNONE && (bDomDec || bMolPBC) && vsite->n_intercg_vsite)
501 /* This is wasting some CPU time as we now do this multiple times
502 * per MD step. But how often do we have vsites with full pbc?
504 pbc_null = set_pbc_dd(&pbc, ePBC, cr != NULL ? cr->dd : NULL, FALSE, box);
513 dd_move_x_vsites(cr->dd, box, x);
516 if (vsite->nthreads == 1)
518 construct_vsites_thread(vsite,
525 #pragma omp parallel num_threads(vsite->nthreads)
527 construct_vsites_thread(vsite,
529 ip, vsite->tdata[gmx_omp_get_thread_num()].ilist,
532 /* Now we can construct the vsites that might depend on other vsites */
533 construct_vsites_thread(vsite,
535 ip, vsite->tdata[vsite->nthreads].ilist,
540 static void spread_vsite2(t_iatom ia[], real a,
541 rvec x[], rvec f[], rvec fshift[],
542 t_pbc *pbc, t_graph *g)
554 svmul(1-a, f[av], fi);
555 svmul( a, f[av], fj);
564 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
566 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
571 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
572 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
580 if (fshift && (siv != CENTRAL || sij != CENTRAL))
582 rvec_inc(fshift[siv], f[av]);
583 rvec_dec(fshift[CENTRAL], fi);
584 rvec_dec(fshift[sij], fj);
587 /* TOTAL: 13 flops */
590 void construct_vsites_mtop(gmx_vsite_t *vsite,
591 gmx_mtop_t *mtop, rvec x[])
594 gmx_molblock_t *molb;
598 for (mb = 0; mb < mtop->nmolblock; mb++)
600 molb = &mtop->molblock[mb];
601 molt = &mtop->moltype[molb->type];
602 for (mol = 0; mol < molb->nmol; mol++)
604 construct_vsites(vsite, x+as, 0.0, NULL,
605 mtop->ffparams.iparams, molt->ilist,
606 epbcNONE, TRUE, NULL, NULL);
607 as += molt->atoms.nr;
612 static void spread_vsite3(t_iatom ia[], real a, real b,
613 rvec x[], rvec f[], rvec fshift[],
614 t_pbc *pbc, t_graph *g)
617 atom_id av, ai, aj, ak;
626 svmul(1-a-b, f[av], fi);
627 svmul( a, f[av], fj);
628 svmul( b, f[av], fk);
638 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ia[1]), di);
640 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), di);
642 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), di);
647 siv = pbc_dx_aiuc(pbc, x[ai], x[av], dx);
648 sij = pbc_dx_aiuc(pbc, x[ai], x[aj], dx);
649 sik = pbc_dx_aiuc(pbc, x[ai], x[ak], dx);
658 if (fshift && (siv != CENTRAL || sij != CENTRAL || sik != CENTRAL))
660 rvec_inc(fshift[siv], f[av]);
661 rvec_dec(fshift[CENTRAL], fi);
662 rvec_dec(fshift[sij], fj);
663 rvec_dec(fshift[sik], fk);
666 /* TOTAL: 20 flops */
669 static void spread_vsite3FD(t_iatom ia[], real a, real b,
670 rvec x[], rvec f[], rvec fshift[],
671 gmx_bool VirCorr, matrix dxdf,
672 t_pbc *pbc, t_graph *g)
674 real fx, fy, fz, c, invl, fproj, a1;
675 rvec xvi, xij, xjk, xix, fv, temp;
676 t_iatom av, ai, aj, ak;
677 int svi, sji, skj, d;
684 copy_rvec(f[av], fv);
686 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
687 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
690 /* xix goes from i to point x on the line jk */
691 xix[XX] = xij[XX]+a*xjk[XX];
692 xix[YY] = xij[YY]+a*xjk[YY];
693 xix[ZZ] = xij[ZZ]+a*xjk[ZZ];
696 invl = gmx_invsqrt(iprod(xix, xix));
700 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
702 temp[XX] = c*(fv[XX]-fproj*xix[XX]);
703 temp[YY] = c*(fv[YY]-fproj*xix[YY]);
704 temp[ZZ] = c*(fv[ZZ]-fproj*xix[ZZ]);
707 /* c is already calculated in constr_vsite3FD
708 storing c somewhere will save 26 flops! */
711 f[ai][XX] += fv[XX] - temp[XX];
712 f[ai][YY] += fv[YY] - temp[YY];
713 f[ai][ZZ] += fv[ZZ] - temp[ZZ];
714 f[aj][XX] += a1*temp[XX];
715 f[aj][YY] += a1*temp[YY];
716 f[aj][ZZ] += a1*temp[ZZ];
717 f[ak][XX] += a*temp[XX];
718 f[ak][YY] += a*temp[YY];
719 f[ak][ZZ] += a*temp[ZZ];
724 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
726 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
728 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
733 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
740 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
742 rvec_dec(fshift[svi], fv);
743 fshift[CENTRAL][XX] += fv[XX] - (1 + a)*temp[XX];
744 fshift[CENTRAL][YY] += fv[YY] - (1 + a)*temp[YY];
745 fshift[CENTRAL][ZZ] += fv[ZZ] - (1 + a)*temp[ZZ];
746 fshift[ sji][XX] += temp[XX];
747 fshift[ sji][YY] += temp[YY];
748 fshift[ sji][ZZ] += temp[ZZ];
749 fshift[ skj][XX] += a*temp[XX];
750 fshift[ skj][YY] += a*temp[YY];
751 fshift[ skj][ZZ] += a*temp[ZZ];
756 /* When VirCorr=TRUE, the virial for the current forces is not
757 * calculated from the redistributed forces. This means that
758 * the effect of non-linear virtual site constructions on the virial
759 * needs to be added separately. This contribution can be calculated
760 * in many ways, but the simplest and cheapest way is to use
761 * the first constructing atom ai as a reference position in space:
762 * subtract (xv-xi)*fv and add (xj-xi)*fj + (xk-xi)*fk.
767 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
769 for (i = 0; i < DIM; i++)
771 for (j = 0; j < DIM; j++)
773 /* As xix is a linear combination of j and k, use that here */
774 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
779 /* TOTAL: 61 flops */
782 static void spread_vsite3FAD(t_iatom ia[], real a, real b,
783 rvec x[], rvec f[], rvec fshift[],
784 gmx_bool VirCorr, matrix dxdf,
785 t_pbc *pbc, t_graph *g)
787 rvec xvi, xij, xjk, xperp, Fpij, Fppp, fv, f1, f2, f3;
788 real a1, b1, c1, c2, invdij, invdij2, invdp, fproj;
789 t_iatom av, ai, aj, ak;
790 int svi, sji, skj, d;
797 copy_rvec(f[ia[1]], fv);
799 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
800 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
803 invdij = gmx_invsqrt(iprod(xij, xij));
804 invdij2 = invdij * invdij;
805 c1 = iprod(xij, xjk) * invdij2;
806 xperp[XX] = xjk[XX] - c1*xij[XX];
807 xperp[YY] = xjk[YY] - c1*xij[YY];
808 xperp[ZZ] = xjk[ZZ] - c1*xij[ZZ];
809 /* xperp in plane ijk, perp. to ij */
810 invdp = gmx_invsqrt(iprod(xperp, xperp));
815 /* a1, b1 and c1 are already calculated in constr_vsite3FAD
816 storing them somewhere will save 45 flops! */
818 fproj = iprod(xij, fv)*invdij2;
819 svmul(fproj, xij, Fpij); /* proj. f on xij */
820 svmul(iprod(xperp, fv)*invdp*invdp, xperp, Fppp); /* proj. f on xperp */
821 svmul(b1*fproj, xperp, f3);
824 rvec_sub(fv, Fpij, f1); /* f1 = f - Fpij */
825 rvec_sub(f1, Fppp, f2); /* f2 = f - Fpij - Fppp */
826 for (d = 0; (d < DIM); d++)
834 f[ai][XX] += fv[XX] - f1[XX] + c1*f2[XX] + f3[XX];
835 f[ai][YY] += fv[YY] - f1[YY] + c1*f2[YY] + f3[YY];
836 f[ai][ZZ] += fv[ZZ] - f1[ZZ] + c1*f2[ZZ] + f3[ZZ];
837 f[aj][XX] += f1[XX] - c2*f2[XX] - f3[XX];
838 f[aj][YY] += f1[YY] - c2*f2[YY] - f3[YY];
839 f[aj][ZZ] += f1[ZZ] - c2*f2[ZZ] - f3[ZZ];
847 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
849 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
851 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
856 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
863 if (fshift && (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL))
865 rvec_dec(fshift[svi], fv);
866 fshift[CENTRAL][XX] += fv[XX] - f1[XX] - (1-c1)*f2[XX] + f3[XX];
867 fshift[CENTRAL][YY] += fv[YY] - f1[YY] - (1-c1)*f2[YY] + f3[YY];
868 fshift[CENTRAL][ZZ] += fv[ZZ] - f1[ZZ] - (1-c1)*f2[ZZ] + f3[ZZ];
869 fshift[ sji][XX] += f1[XX] - c1 *f2[XX] - f3[XX];
870 fshift[ sji][YY] += f1[YY] - c1 *f2[YY] - f3[YY];
871 fshift[ sji][ZZ] += f1[ZZ] - c1 *f2[ZZ] - f3[ZZ];
872 fshift[ skj][XX] += f2[XX];
873 fshift[ skj][YY] += f2[YY];
874 fshift[ skj][ZZ] += f2[ZZ];
882 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
884 for (i = 0; i < DIM; i++)
886 for (j = 0; j < DIM; j++)
888 /* Note that xik=xij+xjk, so we have to add xij*f2 */
891 + xij[i]*(f1[j] + (1 - c2)*f2[j] - f3[j])
897 /* TOTAL: 113 flops */
900 static void spread_vsite3OUT(t_iatom ia[], real a, real b, real c,
901 rvec x[], rvec f[], rvec fshift[],
902 gmx_bool VirCorr, matrix dxdf,
903 t_pbc *pbc, t_graph *g)
905 rvec xvi, xij, xik, fv, fj, fk;
907 atom_id av, ai, aj, ak;
916 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
917 ski = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
920 copy_rvec(f[av], fv);
927 fj[XX] = a*fv[XX] - xik[ZZ]*cfy + xik[YY]*cfz;
928 fj[YY] = xik[ZZ]*cfx + a*fv[YY] - xik[XX]*cfz;
929 fj[ZZ] = -xik[YY]*cfx + xik[XX]*cfy + a*fv[ZZ];
931 fk[XX] = b*fv[XX] + xij[ZZ]*cfy - xij[YY]*cfz;
932 fk[YY] = -xij[ZZ]*cfx + b*fv[YY] + xij[XX]*cfz;
933 fk[ZZ] = xij[YY]*cfx - xij[XX]*cfy + b*fv[ZZ];
936 f[ai][XX] += fv[XX] - fj[XX] - fk[XX];
937 f[ai][YY] += fv[YY] - fj[YY] - fk[YY];
938 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
945 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
947 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
949 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
954 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
961 if (fshift && (svi != CENTRAL || sji != CENTRAL || ski != CENTRAL))
963 rvec_dec(fshift[svi], fv);
964 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX];
965 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY];
966 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ];
967 rvec_inc(fshift[sji], fj);
968 rvec_inc(fshift[ski], fk);
976 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
978 for (i = 0; i < DIM; i++)
980 for (j = 0; j < DIM; j++)
982 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j];
987 /* TOTAL: 54 flops */
990 static void spread_vsite4FD(t_iatom ia[], real a, real b, real c,
991 rvec x[], rvec f[], rvec fshift[],
992 gmx_bool VirCorr, matrix dxdf,
993 t_pbc *pbc, t_graph *g)
995 real d, invl, fproj, a1;
996 rvec xvi, xij, xjk, xjl, xix, fv, temp;
997 atom_id av, ai, aj, ak, al;
999 int svi, sji, skj, slj, m;
1007 sji = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1008 skj = pbc_rvec_sub(pbc, x[ak], x[aj], xjk);
1009 slj = pbc_rvec_sub(pbc, x[al], x[aj], xjl);
1012 /* xix goes from i to point x on the plane jkl */
1013 for (m = 0; m < DIM; m++)
1015 xix[m] = xij[m] + a*xjk[m] + b*xjl[m];
1019 invl = gmx_invsqrt(iprod(xix, xix));
1021 /* 4 + ?10? flops */
1023 copy_rvec(f[av], fv);
1025 fproj = iprod(xix, fv)*invl*invl; /* = (xix . f)/(xix . xix) */
1027 for (m = 0; m < DIM; m++)
1029 temp[m] = d*(fv[m] - fproj*xix[m]);
1033 /* c is already calculated in constr_vsite3FD
1034 storing c somewhere will save 35 flops! */
1037 for (m = 0; m < DIM; m++)
1039 f[ai][m] += fv[m] - temp[m];
1040 f[aj][m] += a1*temp[m];
1041 f[ak][m] += a*temp[m];
1042 f[al][m] += b*temp[m];
1048 ivec_sub(SHIFT_IVEC(g, ia[1]), SHIFT_IVEC(g, ai), di);
1050 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1052 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, aj), di);
1054 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, aj), di);
1059 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1067 (svi != CENTRAL || sji != CENTRAL || skj != CENTRAL || slj != CENTRAL))
1069 rvec_dec(fshift[svi], fv);
1070 for (m = 0; m < DIM; m++)
1072 fshift[CENTRAL][m] += fv[m] - (1 + a + b)*temp[m];
1073 fshift[ sji][m] += temp[m];
1074 fshift[ skj][m] += a*temp[m];
1075 fshift[ slj][m] += b*temp[m];
1084 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1086 for (i = 0; i < DIM; i++)
1088 for (j = 0; j < DIM; j++)
1090 dxdf[i][j] += -xiv[i]*fv[j] + xix[i]*temp[j];
1095 /* TOTAL: 77 flops */
1099 static void spread_vsite4FDN(t_iatom ia[], real a, real b, real c,
1100 rvec x[], rvec f[], rvec fshift[],
1101 gmx_bool VirCorr, matrix dxdf,
1102 t_pbc *pbc, t_graph *g)
1104 rvec xvi, xij, xik, xil, ra, rb, rja, rjb, rab, rm, rt;
1105 rvec fv, fj, fk, fl;
1109 int av, ai, aj, ak, al;
1110 int svi, sij, sik, sil;
1112 /* DEBUG: check atom indices */
1119 copy_rvec(f[av], fv);
1121 sij = pbc_rvec_sub(pbc, x[aj], x[ai], xij);
1122 sik = pbc_rvec_sub(pbc, x[ak], x[ai], xik);
1123 sil = pbc_rvec_sub(pbc, x[al], x[ai], xil);
1136 rvec_sub(ra, xij, rja);
1137 rvec_sub(rb, xij, rjb);
1138 rvec_sub(rb, ra, rab);
1141 cprod(rja, rjb, rm);
1144 invrm = gmx_invsqrt(norm2(rm));
1145 denom = invrm*invrm;
1148 cfx = c*invrm*fv[XX];
1149 cfy = c*invrm*fv[YY];
1150 cfz = c*invrm*fv[ZZ];
1161 fj[XX] = ( -rm[XX]*rt[XX]) * cfx + ( rab[ZZ]-rm[YY]*rt[XX]) * cfy + (-rab[YY]-rm[ZZ]*rt[XX]) * cfz;
1162 fj[YY] = (-rab[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( rab[XX]-rm[ZZ]*rt[YY]) * cfz;
1163 fj[ZZ] = ( rab[YY]-rm[XX]*rt[ZZ]) * cfx + (-rab[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1174 fk[XX] = ( -rm[XX]*rt[XX]) * cfx + (-a*rjb[ZZ]-rm[YY]*rt[XX]) * cfy + ( a*rjb[YY]-rm[ZZ]*rt[XX]) * cfz;
1175 fk[YY] = ( a*rjb[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + (-a*rjb[XX]-rm[ZZ]*rt[YY]) * cfz;
1176 fk[ZZ] = (-a*rjb[YY]-rm[XX]*rt[ZZ]) * cfx + ( a*rjb[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1187 fl[XX] = ( -rm[XX]*rt[XX]) * cfx + ( b*rja[ZZ]-rm[YY]*rt[XX]) * cfy + (-b*rja[YY]-rm[ZZ]*rt[XX]) * cfz;
1188 fl[YY] = (-b*rja[ZZ]-rm[XX]*rt[YY]) * cfx + ( -rm[YY]*rt[YY]) * cfy + ( b*rja[XX]-rm[ZZ]*rt[YY]) * cfz;
1189 fl[ZZ] = ( b*rja[YY]-rm[XX]*rt[ZZ]) * cfx + (-b*rja[XX]-rm[YY]*rt[ZZ]) * cfy + ( -rm[ZZ]*rt[ZZ]) * cfz;
1192 f[ai][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1193 f[ai][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1194 f[ai][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1195 rvec_inc(f[aj], fj);
1196 rvec_inc(f[ak], fk);
1197 rvec_inc(f[al], fl);
1202 ivec_sub(SHIFT_IVEC(g, av), SHIFT_IVEC(g, ai), di);
1204 ivec_sub(SHIFT_IVEC(g, aj), SHIFT_IVEC(g, ai), di);
1206 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, ai), di);
1208 ivec_sub(SHIFT_IVEC(g, al), SHIFT_IVEC(g, ai), di);
1213 svi = pbc_rvec_sub(pbc, x[av], x[ai], xvi);
1220 if (fshift && (svi != CENTRAL || sij != CENTRAL || sik != CENTRAL || sil != CENTRAL))
1222 rvec_dec(fshift[svi], fv);
1223 fshift[CENTRAL][XX] += fv[XX] - fj[XX] - fk[XX] - fl[XX];
1224 fshift[CENTRAL][YY] += fv[YY] - fj[YY] - fk[YY] - fl[YY];
1225 fshift[CENTRAL][ZZ] += fv[ZZ] - fj[ZZ] - fk[ZZ] - fl[ZZ];
1226 rvec_inc(fshift[sij], fj);
1227 rvec_inc(fshift[sik], fk);
1228 rvec_inc(fshift[sil], fl);
1236 pbc_rvec_sub(pbc, x[av], x[ai], xiv);
1238 for (i = 0; i < DIM; i++)
1240 for (j = 0; j < DIM; j++)
1242 dxdf[i][j] += -xiv[i]*fv[j] + xij[i]*fj[j] + xik[i]*fk[j] + xil[i]*fl[j];
1247 /* Total: 207 flops (Yuck!) */
1251 static int spread_vsiten(t_iatom ia[], t_iparams ip[],
1252 rvec x[], rvec f[], rvec fshift[],
1253 t_pbc *pbc, t_graph *g)
1261 n3 = 3*ip[ia[0]].vsiten.n;
1263 copy_rvec(x[av], xv);
1265 for (i = 0; i < n3; i += 3)
1270 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, av), di);
1275 siv = pbc_dx_aiuc(pbc, x[ai], xv, dx);
1281 a = ip[ia[i]].vsiten.a;
1282 svmul(a, f[av], fi);
1283 rvec_inc(f[ai], fi);
1284 if (fshift && siv != CENTRAL)
1286 rvec_inc(fshift[siv], fi);
1287 rvec_dec(fshift[CENTRAL], fi);
1296 static int vsite_count(const t_ilist *ilist, int ftype)
1298 if (ftype == F_VSITEN)
1300 return ilist[ftype].nr/3;
1304 return ilist[ftype].nr/(1 + interaction_function[ftype].nratoms);
1308 static void spread_vsite_f_thread(gmx_vsite_t *vsite,
1309 rvec x[], rvec f[], rvec *fshift,
1310 gmx_bool VirCorr, matrix dxdf,
1311 t_iparams ip[], t_ilist ilist[],
1312 t_graph *g, t_pbc *pbc_null)
1316 int i, inc, m, nra, nr, tp, ftype;
1326 bPBCAll = (pbc_null != NULL && !vsite->bHaveChargeGroups);
1328 /* this loop goes backwards to be able to build *
1329 * higher type vsites from lower types */
1332 for (ftype = F_NRE-1; (ftype >= 0); ftype--)
1334 if ((interaction_function[ftype].flags & IF_VSITE) &&
1335 ilist[ftype].nr > 0)
1337 nra = interaction_function[ftype].nratoms;
1339 nr = ilist[ftype].nr;
1340 ia = ilist[ftype].iatoms;
1344 pbc_null2 = pbc_null;
1346 else if (pbc_null != NULL)
1348 vsite_pbc = vsite->vsite_pbc_loc[ftype-F_VSITE2];
1351 for (i = 0; i < nr; )
1353 if (vsite_pbc != NULL)
1355 if (vsite_pbc[i/(1+nra)] > -2)
1357 pbc_null2 = pbc_null;
1367 /* Constants for constructing */
1368 a1 = ip[tp].vsite.a;
1369 /* Construct the vsite depending on type */
1373 spread_vsite2(ia, a1, x, f, fshift, pbc_null2, g);
1376 b1 = ip[tp].vsite.b;
1377 spread_vsite3(ia, a1, b1, x, f, fshift, pbc_null2, g);
1380 b1 = ip[tp].vsite.b;
1381 spread_vsite3FD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1384 b1 = ip[tp].vsite.b;
1385 spread_vsite3FAD(ia, a1, b1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1388 b1 = ip[tp].vsite.b;
1389 c1 = ip[tp].vsite.c;
1390 spread_vsite3OUT(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1393 b1 = ip[tp].vsite.b;
1394 c1 = ip[tp].vsite.c;
1395 spread_vsite4FD(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1398 b1 = ip[tp].vsite.b;
1399 c1 = ip[tp].vsite.c;
1400 spread_vsite4FDN(ia, a1, b1, c1, x, f, fshift, VirCorr, dxdf, pbc_null2, g);
1403 inc = spread_vsiten(ia, ip, x, f, fshift, pbc_null2, g);
1406 gmx_fatal(FARGS, "No such vsite type %d in %s, line %d",
1407 ftype, __FILE__, __LINE__);
1409 clear_rvec(f[ia[1]]);
1411 /* Increment loop variables */
1419 void spread_vsite_f(gmx_vsite_t *vsite,
1420 rvec x[], rvec f[], rvec *fshift,
1421 gmx_bool VirCorr, matrix vir,
1422 t_nrnb *nrnb, t_idef *idef,
1423 int ePBC, gmx_bool bMolPBC, t_graph *g, matrix box,
1426 t_pbc pbc, *pbc_null;
1429 /* We only need to do pbc when we have inter-cg vsites */
1430 if ((DOMAINDECOMP(cr) || bMolPBC) && vsite->n_intercg_vsite)
1432 /* This is wasting some CPU time as we now do this multiple times
1433 * per MD step. But how often do we have vsites with full pbc?
1435 pbc_null = set_pbc_dd(&pbc, ePBC, cr->dd, FALSE, box);
1442 if (DOMAINDECOMP(cr))
1444 dd_clear_f_vsites(cr->dd, f);
1447 if (vsite->nthreads == 1)
1449 spread_vsite_f_thread(vsite,
1451 VirCorr, vsite->tdata[0].dxdf,
1452 idef->iparams, idef->il,
1457 /* First spread the vsites that might depend on other vsites */
1458 spread_vsite_f_thread(vsite,
1460 VirCorr, vsite->tdata[vsite->nthreads].dxdf,
1462 vsite->tdata[vsite->nthreads].ilist,
1465 #pragma omp parallel num_threads(vsite->nthreads)
1470 thread = gmx_omp_get_thread_num();
1472 if (thread == 0 || fshift == NULL)
1480 fshift_t = vsite->tdata[thread].fshift;
1482 for (i = 0; i < SHIFTS; i++)
1484 clear_rvec(fshift_t[i]);
1488 spread_vsite_f_thread(vsite,
1490 VirCorr, vsite->tdata[thread].dxdf,
1492 vsite->tdata[thread].ilist,
1500 for (th = 1; th < vsite->nthreads; th++)
1502 for (i = 0; i < SHIFTS; i++)
1504 rvec_inc(fshift[i], vsite->tdata[th].fshift[i]);
1514 for (th = 0; th < (vsite->nthreads == 1 ? 1 : vsite->nthreads+1); th++)
1516 for (i = 0; i < DIM; i++)
1518 for (j = 0; j < DIM; j++)
1520 vir[i][j] += -0.5*vsite->tdata[th].dxdf[i][j];
1526 if (DOMAINDECOMP(cr))
1528 dd_move_f_vsites(cr->dd, f, fshift);
1531 inc_nrnb(nrnb, eNR_VSITE2, vsite_count(idef->il, F_VSITE2));
1532 inc_nrnb(nrnb, eNR_VSITE3, vsite_count(idef->il, F_VSITE3));
1533 inc_nrnb(nrnb, eNR_VSITE3FD, vsite_count(idef->il, F_VSITE3FD));
1534 inc_nrnb(nrnb, eNR_VSITE3FAD, vsite_count(idef->il, F_VSITE3FAD));
1535 inc_nrnb(nrnb, eNR_VSITE3OUT, vsite_count(idef->il, F_VSITE3OUT));
1536 inc_nrnb(nrnb, eNR_VSITE4FD, vsite_count(idef->il, F_VSITE4FD));
1537 inc_nrnb(nrnb, eNR_VSITE4FDN, vsite_count(idef->il, F_VSITE4FDN));
1538 inc_nrnb(nrnb, eNR_VSITEN, vsite_count(idef->il, F_VSITEN));
1541 static int *atom2cg(t_block *cgs)
1545 snew(a2cg, cgs->index[cgs->nr]);
1546 for (cg = 0; cg < cgs->nr; cg++)
1548 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
1557 static int count_intercg_vsite(gmx_mtop_t *mtop,
1558 gmx_bool *bHaveChargeGroups)
1560 int mb, mt, ftype, nral, i, cg, a;
1561 gmx_molblock_t *molb;
1562 gmx_moltype_t *molt;
1566 int n_intercg_vsite;
1568 *bHaveChargeGroups = FALSE;
1570 n_intercg_vsite = 0;
1571 for (mb = 0; mb < mtop->nmolblock; mb++)
1573 molb = &mtop->molblock[mb];
1574 molt = &mtop->moltype[molb->type];
1576 if (molt->cgs.nr < molt->atoms.nr)
1578 *bHaveChargeGroups = TRUE;
1581 a2cg = atom2cg(&molt->cgs);
1582 for (ftype = 0; ftype < F_NRE; ftype++)
1584 if (interaction_function[ftype].flags & IF_VSITE)
1587 il = &molt->ilist[ftype];
1589 for (i = 0; i < il->nr; i += 1+nral)
1592 for (a = 1; a < nral; a++)
1594 if (a2cg[ia[1+a]] != cg)
1596 n_intercg_vsite += molb->nmol;
1606 return n_intercg_vsite;
1609 static int **get_vsite_pbc(t_iparams *iparams, t_ilist *ilist,
1610 t_atom *atom, t_mdatoms *md,
1611 t_block *cgs, int *a2cg)
1613 int ftype, nral, i, j, vsi, vsite, cg_v, cg_c, a, nc3 = 0;
1616 int **vsite_pbc, *vsite_pbc_f;
1618 gmx_bool bViteOnlyCG_and_FirstAtom;
1620 /* Make an array that tells if the pbc of an atom is set */
1621 snew(pbc_set, cgs->index[cgs->nr]);
1622 /* PBC is set for all non vsites */
1623 for (a = 0; a < cgs->index[cgs->nr]; a++)
1625 if ((atom && atom[a].ptype != eptVSite) ||
1626 (md && md->ptype[a] != eptVSite))
1632 snew(vsite_pbc, F_VSITEN-F_VSITE2+1);
1634 for (ftype = 0; ftype < F_NRE; ftype++)
1636 if (interaction_function[ftype].flags & IF_VSITE)
1642 snew(vsite_pbc[ftype-F_VSITE2], il->nr/(1+nral));
1643 vsite_pbc_f = vsite_pbc[ftype-F_VSITE2];
1651 /* A value of -2 signals that this vsite and its contructing
1652 * atoms are all within the same cg, so no pbc is required.
1654 vsite_pbc_f[vsi] = -2;
1655 /* Check if constructing atoms are outside the vsite's cg */
1657 if (ftype == F_VSITEN)
1659 nc3 = 3*iparams[ia[i]].vsiten.n;
1660 for (j = 0; j < nc3; j += 3)
1662 if (a2cg[ia[i+j+2]] != cg_v)
1664 vsite_pbc_f[vsi] = -1;
1670 for (a = 1; a < nral; a++)
1672 if (a2cg[ia[i+1+a]] != cg_v)
1674 vsite_pbc_f[vsi] = -1;
1678 if (vsite_pbc_f[vsi] == -1)
1680 /* Check if this is the first processed atom of a vsite only cg */
1681 bViteOnlyCG_and_FirstAtom = TRUE;
1682 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1684 /* Non-vsites already have pbc set, so simply check for pbc_set */
1687 bViteOnlyCG_and_FirstAtom = FALSE;
1691 if (bViteOnlyCG_and_FirstAtom)
1693 /* First processed atom of a vsite only charge group.
1694 * The pbc of the input coordinates to construct_vsites
1695 * should be preserved.
1697 vsite_pbc_f[vsi] = vsite;
1699 else if (cg_v != a2cg[ia[1+i+1]])
1701 /* This vsite has a different charge group index
1702 * than it's first constructing atom
1703 * and the charge group has more than one atom,
1704 * search for the first normal particle
1705 * or vsite that already had its pbc defined.
1706 * If nothing is found, use full pbc for this vsite.
1708 for (a = cgs->index[cg_v]; a < cgs->index[cg_v+1]; a++)
1710 if (a != vsite && pbc_set[a])
1712 vsite_pbc_f[vsi] = a;
1715 fprintf(debug, "vsite %d match pbc with atom %d\n",
1723 fprintf(debug, "vsite atom %d cg %d - %d pbc atom %d\n",
1724 vsite+1, cgs->index[cg_v]+1, cgs->index[cg_v+1],
1725 vsite_pbc_f[vsi]+1);
1729 if (ftype == F_VSITEN)
1731 /* The other entries in vsite_pbc_f are not used for center vsites */
1739 /* This vsite now has its pbc defined */
1751 gmx_vsite_t *init_vsite(gmx_mtop_t *mtop, t_commrec *cr,
1752 gmx_bool bSerial_NoPBC)
1758 gmx_moltype_t *molt;
1761 /* check if there are vsites */
1763 for (i = 0; i < F_NRE; i++)
1765 if (interaction_function[i].flags & IF_VSITE)
1767 nvsite += gmx_mtop_ftype_count(mtop, i);
1778 vsite->n_intercg_vsite = count_intercg_vsite(mtop,
1779 &vsite->bHaveChargeGroups);
1781 /* If we don't have charge groups, the vsite follows its own pbc */
1782 if (!bSerial_NoPBC &&
1783 vsite->bHaveChargeGroups &&
1784 vsite->n_intercg_vsite > 0 && DOMAINDECOMP(cr))
1786 vsite->nvsite_pbc_molt = mtop->nmoltype;
1787 snew(vsite->vsite_pbc_molt, vsite->nvsite_pbc_molt);
1788 for (mt = 0; mt < mtop->nmoltype; mt++)
1790 molt = &mtop->moltype[mt];
1791 /* Make an atom to charge group index */
1792 a2cg = atom2cg(&molt->cgs);
1793 vsite->vsite_pbc_molt[mt] = get_vsite_pbc(mtop->ffparams.iparams,
1795 molt->atoms.atom, NULL,
1800 snew(vsite->vsite_pbc_loc_nalloc, F_VSITEN-F_VSITE2+1);
1801 snew(vsite->vsite_pbc_loc, F_VSITEN-F_VSITE2+1);
1806 vsite->nthreads = 1;
1810 vsite->nthreads = gmx_omp_nthreads_get(emntVSITE);
1814 /* We need one extra thread data structure for the overlap vsites */
1815 snew(vsite->tdata, vsite->nthreads+1);
1818 vsite->th_ind = NULL;
1819 vsite->th_ind_nalloc = 0;
1824 static void prepare_vsite_thread(const t_ilist *ilist,
1825 gmx_vsite_thread_t *vsite_th)
1829 for (ftype = 0; ftype < F_NRE; ftype++)
1831 if (interaction_function[ftype].flags & IF_VSITE)
1833 if (ilist[ftype].nr > vsite_th->ilist[ftype].nalloc)
1835 vsite_th->ilist[ftype].nalloc = over_alloc_large(ilist[ftype].nr);
1836 srenew(vsite_th->ilist[ftype].iatoms, vsite_th->ilist[ftype].nalloc);
1839 vsite_th->ilist[ftype].nr = 0;
1844 void split_vsites_over_threads(const t_ilist *ilist,
1845 const t_mdatoms *mdatoms,
1846 gmx_bool bLimitRange,
1850 int vsite_atom_range, natperthread;
1855 int nral1, inc, i, j;
1857 if (vsite->nthreads == 1)
1863 #pragma omp parallel for num_threads(vsite->nthreads) schedule(static)
1864 for (th = 0; th < vsite->nthreads; th++)
1866 prepare_vsite_thread(ilist, &vsite->tdata[th]);
1868 /* Master threads does the (potential) overlap vsites */
1869 prepare_vsite_thread(ilist, &vsite->tdata[vsite->nthreads]);
1871 /* The current way of distributing the vsites over threads in primitive.
1872 * We divide the atom range 0 - natoms_in_vsite uniformly over threads,
1873 * without taking into account how the vsites are distributed.
1874 * Without domain decomposition we bLimitRange=TRUE and we at least
1875 * tighten the upper bound of the range (useful for common systems
1876 * such as a vsite-protein in 3-site water).
1880 vsite_atom_range = -1;
1881 for (ftype = 0; ftype < F_NRE; ftype++)
1883 if ((interaction_function[ftype].flags & IF_VSITE) &&
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]);
1901 vsite_atom_range = mdatoms->homenr;
1903 natperthread = (vsite_atom_range + vsite->nthreads - 1)/vsite->nthreads;
1907 fprintf(debug, "virtual site thread dist: natoms %d, range %d, natperthread %d\n", mdatoms->nr, vsite_atom_range, natperthread);
1910 /* To simplify the vsite assignment, we make an index which tells us
1911 * to which thread particles, both non-vsites and vsites, are assigned.
1913 if (mdatoms->nr > vsite->th_ind_nalloc)
1915 vsite->th_ind_nalloc = over_alloc_large(mdatoms->nr);
1916 srenew(vsite->th_ind, vsite->th_ind_nalloc);
1918 th_ind = vsite->th_ind;
1920 for (i = 0; i < mdatoms->nr; i++)
1922 if (mdatoms->ptype[i] == eptVSite)
1924 /* vsites are not assigned to a thread yet */
1929 /* assign non-vsite particles to thread th */
1932 if (i == (th + 1)*natperthread && th < vsite->nthreads)
1938 for (ftype = 0; ftype < F_NRE; ftype++)
1940 if ((interaction_function[ftype].flags & IF_VSITE) &&
1943 nral1 = 1 + NRAL(ftype);
1945 iat = ilist[ftype].iatoms;
1946 for (i = 0; i < ilist[ftype].nr; )
1948 th = iat[1+i]/natperthread;
1949 /* We would like to assign this vsite the thread th,
1950 * but it might depend on atoms outside the atom range of th
1951 * or on another vsite not assigned to thread th.
1953 if (ftype != F_VSITEN)
1955 for (j = i+2; j < i+nral1; j++)
1957 if (th_ind[iat[j]] != th)
1959 /* Some constructing atoms are not assigned to
1960 * thread th, move this vsite to a separate batch.
1962 th = vsite->nthreads;
1969 for (j = i+2; j < i+inc; j += 3)
1971 if (th_ind[iat[j]] != th)
1973 th = vsite->nthreads;
1977 /* Copy this vsite to the thread data struct of thread th */
1978 il_th = &vsite->tdata[th].ilist[ftype];
1979 for (j = i; j < i+inc; j++)
1981 il_th->iatoms[il_th->nr++] = iat[j];
1983 /* Update this vsite's thread index entry */
1984 th_ind[iat[1+i]] = th;
1993 for (ftype = 0; ftype < F_NRE; ftype++)
1995 if ((interaction_function[ftype].flags & IF_VSITE) &&
1996 ilist[ftype].nr > 0)
1998 fprintf(debug, "%-20s thread dist:",
1999 interaction_function[ftype].longname);
2000 for (th = 0; th < vsite->nthreads+1; th++)
2002 fprintf(debug, " %4d", vsite->tdata[th].ilist[ftype].nr);
2004 fprintf(debug, "\n");
2010 void set_vsite_top(gmx_vsite_t *vsite, gmx_localtop_t *top, t_mdatoms *md,
2015 if (vsite->n_intercg_vsite > 0)
2017 if (vsite->bHaveChargeGroups)
2019 /* Make an atom to charge group index */
2020 a2cg = atom2cg(&top->cgs);
2021 vsite->vsite_pbc_loc = get_vsite_pbc(top->idef.iparams,
2022 top->idef.il, NULL, md,
2029 if (vsite->nthreads > 1)
2031 if (vsite->bHaveChargeGroups)
2033 gmx_fatal(FARGS, "The combination of threading, virtual sites and charge groups is not implemented");
2036 split_vsites_over_threads(top->idef.il, md, !DOMAINDECOMP(cr), vsite);